1 #include "heatsystem.h" 3 #include "libmesh/dirichlet_boundaries.h" 4 #include "libmesh/dof_map.h" 5 #include "libmesh/fe_base.h" 6 #include "libmesh/fe_interface.h" 7 #include "libmesh/fem_context.h" 8 #include "libmesh/getpot.h" 9 #include "libmesh/mesh.h" 10 #include "libmesh/quadrature.h" 11 #include "libmesh/string_to_enum.h" 12 #include "libmesh/zero_function.h" 13 #include "libmesh/elem.h" 19 T_var = this->add_variable(
"T", static_cast<Order>(_fe_order),
20 Utility::string_to_enum<FEFamily>(_fe_family));
22 #ifdef LIBMESH_ENABLE_DIRICHLET 23 const unsigned int dim = this->get_mesh().mesh_dimension();
27 std::set<boundary_id_type> nonyplus_bdys(all_ids, all_ids+(
dim*2));
29 nonyplus_bdys.erase(yplus_id);
35 this->get_dof_map().add_dirichlet_boundary
37 #endif // LIBMESH_ENABLE_DIRICHLET 43 this->time_evolving(T_var, 1);
50 FEMContext & c = cast_ref<FEMContext &>(context);
52 const std::set<unsigned char> & elem_dims =
55 for (std::set<unsigned char>::const_iterator dim_it =
56 elem_dims.begin(); dim_it != elem_dims.end(); ++dim_it)
58 const unsigned char dim = *dim_it;
69 FEBase * side_fe =
nullptr;
81 FEMContext & c = cast_ref<FEMContext &>(context);
83 const unsigned int mesh_dim =
93 const std::vector<Real> & JxW = fe->
get_JxW();
95 const std::vector<Point> & xyz = fe->
get_xyz();
97 const std::vector<std::vector<Real>> & phi = fe->
get_phi();
99 const std::vector<std::vector<RealGradient>> & dphi = fe->
get_dphi();
116 for (
unsigned int qp=0; qp != n_qpoints; qp++)
123 const Point & p = xyz[qp];
126 const Number u_exact = (mesh_dim == 2) ?
132 const Number forcing = (
dim == mesh_dim) ?
135 const Number JxWxNK = JxW[qp] * -k;
137 for (
unsigned int i=0; i != n_T_dofs; i++)
138 F(i) += JxWxNK * (grad_T * dphi[i][qp] + forcing * phi[i][qp]);
139 if (request_jacobian)
141 const Number JxWxNKxD = JxWxNK *
144 for (
unsigned int i=0; i != n_T_dofs; i++)
145 for (
unsigned int j=0; j != n_T_dofs; ++j)
146 K(i,j) += JxWxNKxD * (dphi[i][qp] * dphi[j][qp]);
150 return request_jacobian;
Real get_elem_solution_derivative() const
The derivative of the current elem_solution w.r.t.
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
This class provides all data required for a physics package (e.g.
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...
ConstFunction that simply returns 0.
const Elem & get_elem() const
Accessor for current Elem object.
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
virtual void init_context(DiffContext &) override
unsigned int n_dof_indices() const
Total number of dof indices on the element.
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
The libMesh namespace provides an interface to certain functionality in the library.
const MeshBase & get_mesh() const
Defines a dense subvector for use in finite element computations.
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Gradient interior_gradient(unsigned int var, unsigned int qp) const
const System & get_system() const
Accessor for associated system.
virtual void init_context(DiffContext &context)
Defines a dense submatrix for use in Finite Element-type computations.
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
virtual_for_inffe const std::vector< Real > & get_JxW() const
This class provides all data required for a physics package (e.g.
unsigned int n_points() const
virtual_for_inffe const std::vector< Point > & get_xyz() const
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
virtual unsigned short dim() const =0
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
unsigned int mesh_dimension() const
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...
A Point defines a location in LIBMESH_DIM dimensional Real space.
This class forms the foundation from which generic finite elements may be derived.
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
const std::set< unsigned char > & elem_dimensions() const
const std::vector< std::vector< OutputShape > > & get_phi() const