Go to the documentation of this file. 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);
31 std::vector<unsigned int> T_only(1, T_var);
36 this->get_dof_map().add_dirichlet_boundary
38 #endif // LIBMESH_ENABLE_DIRICHLET
44 this->time_evolving(T_var, 1);
51 FEMContext & c = cast_ref<FEMContext &>(context);
53 const std::set<unsigned char> & elem_dims =
56 for (std::set<unsigned char>::const_iterator dim_it =
57 elem_dims.begin(); dim_it != elem_dims.end(); ++dim_it)
59 const unsigned char dim = *dim_it;
78 FEMContext & c = cast_ref<FEMContext &>(context);
80 const unsigned int mesh_dim =
90 const std::vector<Real> & JxW = fe->
get_JxW();
92 const std::vector<Point> & xyz = fe->
get_xyz();
94 const std::vector<std::vector<Real>> & phi = fe->
get_phi();
96 const std::vector<std::vector<RealGradient>> & dphi = fe->
get_dphi();
113 for (
unsigned int qp=0; qp != n_qpoints; qp++)
120 const Point & p = xyz[qp];
123 const Number u_exact = (mesh_dim == 2) ?
129 const Number forcing = (
dim == mesh_dim) ?
132 const Number JxWxNK = JxW[qp] * -k;
134 for (
unsigned int i=0; i != n_T_dofs; i++)
135 F(i) += JxWxNK * (grad_T * dphi[i][qp] + forcing * phi[i][qp]);
136 if (request_jacobian)
138 const Number JxWxNKxD = JxWxNK *
141 for (
unsigned int i=0; i != n_T_dofs; i++)
142 for (
unsigned int j=0; j != n_T_dofs; ++j)
143 K(i,j) += JxWxNKxD * (dphi[i][qp] * dphi[j][qp]);
147 return request_jacobian;
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.
virtual void init_context(DiffContext &) override
const std::vector< Point > & get_xyz() const
virtual void init_context(DiffContext &context)
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.
virtual unsigned short dim() const =0
Real get_elem_solution_derivative() const
The derivative of the current elem_solution w.r.t.
This class forms the foundation from which generic finite elements may be derived.
const System & get_system() const
Accessor for associated system.
const std::vector< Real > & get_JxW() const
const Elem & get_elem() const
Accessor for current Elem object.
unsigned int mesh_dimension() const
const std::vector< std::vector< OutputGradient > > & get_dphi() const
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
const MeshBase & get_mesh() const
A Point defines a location in LIBMESH_DIM dimensional Real space.
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
unsigned int n_dof_indices() const
Total number of dof indices on the element.
ConstFunction that simply returns 0.
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
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
const std::set< unsigned char > & elem_dimensions() const
const std::vector< std::vector< OutputShape > > & get_phi() const
This class provides all data required for a physics package (e.g.