Go to the documentation of this file.
4 #include "libmesh/getpot.h"
5 #include "libmesh/fe_base.h"
6 #include "libmesh/quadrature.h"
7 #include "libmesh/string_to_enum.h"
8 #include "libmesh/parallel.h"
9 #include "libmesh/fem_context.h"
20 this->add_variable (
"T", static_cast<Order>(_fe_order),
21 Utility::string_to_enum<FEFamily>(_fe_family));
23 GetPot infile(
"l-shaped.in");
24 parameters.push_back(infile(
"alpha_1", 1.0));
25 parameters.push_back(infile(
"alpha_2", 1.0));
31 this->time_evolving(T_var, 1);
36 FEMContext & c = cast_ref<FEMContext &>(context);
40 FEBase * elem_fe =
nullptr;
45 FEBase * side_fe =
nullptr;
52 #define optassert(X) {if (!(X)) libmesh_error_msg("Assertion " #X " failed.");}
61 FEMContext & c = cast_ref<FEMContext &>(context);
65 FEBase * elem_fe =
nullptr;
69 const std::vector<Real> & JxW = elem_fe->
get_JxW();
72 const std::vector<std::vector<RealGradient>> & dphi = elem_fe->
get_dphi();
89 for (
unsigned int qp=0; qp != n_qpoints; qp++)
95 for (
unsigned int i=0; i != n_T_dofs; i++)
96 F(i) += JxW[qp] * (parameters[0] + (2.*parameters[1])) * (grad_T * dphi[i][qp]);
98 for (
unsigned int i=0; i != n_T_dofs; i++)
99 for (
unsigned int j=0; j != n_T_dofs; ++j)
101 K(i,j) += JxW[qp] * (parameters[0] + (2.*parameters[1])) * (dphi[i][qp] * dphi[j][qp]);
114 FEMContext & c = cast_ref<FEMContext &>(context);
118 FEBase * side_fe =
nullptr;
122 const std::vector<Real> & JxW = side_fe->
get_JxW();
125 const std::vector<std::vector<Real>> & phi = side_fe->
get_phi();
128 const std::vector<Point > & qside_point = side_fe->
get_xyz();
141 for (
unsigned int qp=0; qp != n_qpoints; qp++)
150 for (
unsigned int i=0; i != n_T_dofs; i++)
151 F(i) += JxW[qp] * penalty * (T - u_dirichlet) * phi[i][qp];
153 for (
unsigned int i=0; i != n_T_dofs; i++)
154 for (
unsigned int j=0; j != n_T_dofs; ++j)
156 K(i,j) += JxW[qp] * penalty * phi[i][qp] * phi[j][qp];
167 const Real x1 = p(0);
168 const Real x2 = p(1);
170 Real theta = atan2(x2, x1);
176 return (parameters[0] + (2.*parameters[1])) *
pow(x1*x1 + x2*x2, 1./3.)*sin(2./3.*theta);
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
const std::vector< Point > & get_xyz() const
Number(* exact_solution)(const Point &p, const Parameters &, const std::string &, const std::string &)
void compute_jacobian(const NumericVector< Number > &, SparseMatrix< Number > &J, NonlinearImplicitSystem &system)
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
Defines a dense submatrix for use in Finite Element-type computations.
unsigned int n_points() const
The libMesh namespace provides an interface to certain functionality in the library.
This class forms the foundation from which generic finite elements may be derived.
const std::vector< Real > & get_JxW() const
static const Real TOLERANCE
const std::vector< std::vector< OutputGradient > > & get_dphi() const
LaplaceExactSolution exact_solution
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
A Point defines a location in LIBMESH_DIM dimensional Real space.
unsigned int n_dof_indices() const
Total number of dof indices on the element.
double pow(double a, int b)
This class provides all data required for a physics package (e.g.
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh.
Gradient interior_gradient(unsigned int var, unsigned int qp) const
virtual void init_context(DiffContext &context)
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Number side_value(unsigned int var, unsigned int qp) const
const QBase & get_side_qrule() const
Accessor for element side quadrature rule for the dimension of the current _elem.
virtual bool side_constraint(bool request_jacobian, DiffContext &context)
Adds the constraint contribution on side of elem to elem_residual.
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in th...
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
const std::vector< std::vector< OutputShape > > & get_phi() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
This class provides all data required for a physics package (e.g.