Go to the documentation of this file.
   20 #ifndef LIBMESH_FEM_CONTEXT_H 
   21 #define LIBMESH_FEM_CONTEXT_H 
   24 #include "libmesh/diff_context.h" 
   25 #include "libmesh/id_types.h" 
   26 #include "libmesh/fe_type.h" 
   27 #include "libmesh/fe_base.h" 
   28 #include "libmesh/vector_value.h" 
   30 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
   31 #include "libmesh/tensor_value.h" 
   44 template <
typename T> 
class FEGenericBase;
 
   45 typedef FEGenericBase<Real> 
FEBase;
 
   48 template <
typename T> 
class NumericVector;
 
  107 #ifdef LIBMESH_ENABLE_DEPRECATED 
  164 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  189 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 
  239 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  264 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 
  274   template<
typename OutputShape>
 
  276   { this->get_element_fe<OutputShape>(var,fe,this->
get_dim()); }
 
  293   template<
typename OutputShape>
 
  295                        unsigned short dim ) 
const;
 
  311   template<
typename OutputShape>
 
  313   { this->get_side_fe<OutputShape>(var,fe,this->
get_dim()); }
 
  330   template<
typename OutputShape>
 
  332                     unsigned short dim ) 
const;
 
  343   template<
typename OutputShape>
 
  357   template<
typename OutputType>
 
  360                       OutputType & u) 
const;
 
  366   template<
typename OutputType>
 
  369                        std::vector<OutputType> & interior_values_vector) 
const;
 
  377   template<
typename OutputType>
 
  380                   OutputType & u) 
const;
 
  386   template<
typename OutputType>
 
  389                    std::vector<OutputType> & side_values_vector) 
const;
 
  400   template<
typename OutputType>
 
  412   template<
typename OutputType>
 
  414                          OutputType & du) 
const;
 
  422   template<
typename OutputType>
 
  425                           std::vector<OutputType> & interior_gradients_vector) 
const;
 
  433   template<
typename OutputType>
 
  436                      OutputType & du) 
const;
 
  444   template<
typename OutputType>
 
  447                       std::vector<OutputType> & side_gradients_vector) 
const;
 
  458   template<
typename OutputType>
 
  464 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  471   template<
typename OutputType>
 
  474                         OutputType & d2u) 
const;
 
  481   template<
typename OutputType>
 
  484                          std::vector<OutputType> & d2u_vals) 
const;
 
  492   template<
typename OutputType>
 
  495                     OutputType & d2u) 
const;
 
  502   template<
typename OutputType>
 
  505                      std::vector<OutputType> & d2u_vals) 
const;
 
  516   template<
typename OutputType>
 
  522 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 
  529   template<
typename OutputType>
 
  532                      OutputType & u) 
const;
 
  540   template<
typename OutputType>
 
  543                               OutputType & u) 
const;
 
  550   template<
typename OutputType>
 
  553                  OutputType & u) 
const;
 
  559   template<
typename OutputType>
 
  562                   OutputType & u) 
const;
 
  569   template<
typename OutputType>
 
  572                       OutputType & u) 
const;
 
  578   template<
typename OutputType>
 
  581                   OutputType & u) 
const;
 
  587   template<
typename OutputType>
 
  590                    OutputType & u) 
const;
 
  598   template<
typename OutputType>
 
  601                             OutputType & u) 
const;
 
  609   template<
typename OutputType>
 
  612                         OutputType & u) 
const;
 
  623   template<
typename OutputType>
 
  635   template<
typename OutputType>
 
  638                                OutputType & grad_u) 
const;
 
  646   template<
typename OutputType>
 
  649                            OutputType & grad_u) 
const;
 
  660   template<
typename OutputType>
 
  666 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  673   template<
typename OutputType>
 
  676                               OutputType & hess_u) 
const;
 
  684   template<
typename OutputType>
 
  687                           OutputType & hess_u) 
const;
 
  698   template<
typename OutputType>
 
  704 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 
  710   template<
typename OutputType>
 
  713                      OutputType & curl_u) 
const;
 
  722   template<
typename OutputType>
 
  732   template<
typename OutputType>
 
  735                     OutputType & div_u) 
const;
 
  774   virtual void elem_fe_reinit(
const std::vector<Point> * 
const pts = 
nullptr);
 
  891   { 
return (this->
_elem != 
nullptr); }
 
  898     return *(this->
_elem); }
 
  905     return *(const_cast<Elem *>(this->
_elem)); }
 
  925   { 
return this->
_dim; }
 
 1030   template<
typename OutputShape>
 
 1050 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 
 1055   template<
typename OutputShape>
 
 1057                                           const FEType fe_type ) 
const;
 
 1077   template <
typename OutputType>
 
 1106   template<
typename OutputType,
 
 1109   void some_value(
unsigned int var, 
unsigned int qp, OutputType & u) 
const;
 
 1115   template<
typename OutputType,
 
 1118   void some_gradient(
unsigned int var, 
unsigned int qp, OutputType & u) 
const;
 
 1120 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
 1125   template<
typename OutputType,
 
 1128   void some_hessian(
unsigned int var, 
unsigned int qp, OutputType & u) 
const;
 
 1136   std::vector<std::map<FEType, std::unique_ptr<FEAbstract>>> 
_element_fe;
 
 1137   std::vector<std::map<FEType, std::unique_ptr<FEAbstract>>> 
_side_fe;
 
 1138   std::map<FEType, std::unique_ptr<FEAbstract>> 
_edge_fe;
 
 1244 template<
typename OutputShape>
 
 1247                                  unsigned short dim )
 const 
 1262 template<
typename OutputShape>
 
 1265                               unsigned short dim )
 const 
 1269   fe = cast_ptr<FEGenericBase<OutputShape> *>( (
_side_fe_var[
dim][var] ) );
 
 1280 template<
typename OutputShape>
 
 1285   fe = cast_ptr<FEGenericBase<OutputShape> *>( 
_edge_fe_var[var] );
 
 1298 #endif // LIBMESH_FEM_CONTEXT_H 
  
The QBase class provides the basic functionality from which various quadrature rules can be derived.
 
Manages consistently variables, degrees of freedom, and coefficient vectors.
 
TensorTools::MakeReal< Rank2Decrement >::type hess_shape
 
Helper nested class for C++03-compatible "template typedef".
 
void init_internal_data(const System &sys)
Helper function used in constructors to set up internal data.
 
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
 
Gradient fixed_point_gradient(unsigned int var, const Point &p) const
 
DiffContext(const System &)
Constructor.
 
The BoundaryInfo class contains information relevant to boundary conditions including storing faces,...
 
Gradient side_gradient(unsigned int var, unsigned int qp) const
 
unsigned char side
Current side for side_* to examine.
 
void elem_position_get()
Uses the geometry of elem to set the coordinate data specified by mesh_*_position configuration.
 
FEGenericBase< grad_shape > grad_base
 
const QBase & get_element_qrule(unsigned short dim) const
Accessor for element interior quadrature rule.
 
std::unique_ptr< FEGenericBase< RealGradient > > _real_grad_fe
 
void(FEMContext::* hess_getter)(unsigned int, hess_base *&, unsigned short) const
 
virtual ~FEMContext()
Destructor.
 
virtual void nonlocal_reinit(Real theta) override
Gives derived classes the opportunity to reinitialize data needed for nonlocal calculations at a new ...
 
FEMContext(const System &sys)
Constructor.
 
std::vector< std::unique_ptr< QBase > > _element_qrule
Quadrature rule for element interior.
 
TensorTools::MakeReal< Rank1Decrement >::type grad_shape
 
void side_accel(unsigned int var, unsigned int qp, OutputType &u) const
 
unsigned int get_mesh_y_var() const
Accessor for y-variable of moving mesh System.
 
virtual void elem_side_reinit(Real theta) override
Resets the current time in the context.
 
std::unique_ptr< QBase > _edge_qrule
Quadrature rules for element edges.
 
std::vector< boundary_id_type > side_boundary_ids() const
Lists the boundary ids found on the current side.
 
void side_gradients(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_gradients_vector) const
Fills a vector with the gradient of the solution variable var at all the quadrature points on the cur...
 
Gradient fixed_side_gradient(unsigned int var, unsigned int qp) const
 
virtual void side_fe_reinit()
Reinitializes side FE objects on the current geometric element.
 
TensorTools::DecrementRank< Rank1Decrement >::type Rank2Decrement
 
FEGenericBase< OutputShape > * build_new_fe(const FEGenericBase< OutputShape > *fe, const Point &p, const Real tolerance=TOLERANCE) const
Helper function to reduce some code duplication in the *_point_* methods.
 
void set_elem(const Elem *e)
Helper function to promote accessor usage.
 
void interior_gradients(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_gradients_vector) const
Fills a vector with the gradient of the solution variable var at all the quadrature points in the cur...
 
The libMesh namespace provides an interface to certain functionality in the library.
 
unsigned char get_side() const
Accessor for current side of Elem object.
 
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
 
This class forms the foundation from which generic finite elements may be derived.
 
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _side_fe
 
FEBase * get_element_fe(unsigned int var) const
Accessor for interior finite element object for scalar-valued variable var for the largest dimension ...
 
unsigned char get_elem_dim() const
 
const Elem & get_elem() const
Accessor for current Elem object.
 
static const Real TOLERANCE
 
void get_edge_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge (3D only!) finite element object for variable var.
 
System * get_mesh_system()
Accessor for moving mesh System.
 
void interior_accel(unsigned int var, unsigned int qp, OutputType &u) const
 
void interior_curl(unsigned int var, unsigned int qp, OutputType &curl_u) const
 
void(FEMContext::* grad_getter)(unsigned int, grad_base *&, unsigned short) const
 
TensorTools::DecrementRank< OutputType >::type Rank1Decrement
 
const BoundaryInfo & _boundary_info
Saved reference to BoundaryInfo on the mesh for this System.
 
Number point_value(unsigned int var, const Point &p) const
 
void set_custom_solution(const NumericVector< Number > *custom_sol)
Set a NumericVector to be used in place of current_local_solution for calculating elem_solution.
 
int _extra_quadrature_order
The extra quadrature order for this context.
 
Number fixed_point_value(unsigned int var, const Point &p) const
 
unsigned char _elem_dim
Cached dimension of this->_elem.
 
void some_gradient(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_gradient methods.
 
std::vector< std::vector< FEAbstract * > > _side_fe_var
 
void use_unweighted_quadrature_rules(int extra_quadrature_order=0)
Use quadrature rules designed to exactly integrate unweighted undistorted basis functions,...
 
Number fixed_side_value(unsigned int var, unsigned int qp) const
 
void point_rate(unsigned int var, const Point &p, OutputType &u) const
 
FEGenericBase< OutputShape > * cached_fe(const unsigned int elem_dim, const FEType fe_type) const
 
void point_curl(unsigned int var, const Point &p, OutputType &curl_u, const Real tolerance=TOLERANCE) const
 
void interior_rate(unsigned int var, unsigned int qp, OutputType &u) const
 
std::set< unsigned char > _elem_dims
Cached dimensions of elements in the mesh, plus dimension 0 if SCALAR variables are in use.
 
void some_hessian(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_hessian methods.
 
void side_rate(unsigned int var, unsigned int qp, OutputType &u) const
 
const Elem * _elem
Current element for element_* to examine.
 
void set_jacobian_tolerance(Real tol)
Calls set_jacobian_tolerance() on all the FE objects controlled by this class.
 
System * _mesh_sys
System from which to acquire moving mesh information.
 
bool has_elem() const
Test for current Elem object.
 
FEGenericBase< Real > FEBase
 
const QBase & get_side_qrule(unsigned short dim) const
Accessor for element side quadrature rule.
 
Tensor interior_hessian(unsigned int var, unsigned int qp) const
 
std::unique_ptr< FEGenericBase< Real > > _real_fe
 
const System * get_mesh_system() const
Accessor for moving mesh System.
 
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _element_fe
Finite element objects for each variable's interior, sides and edges.
 
Gradient point_gradient(unsigned int var, const Point &p) const
 
FEBase * get_side_fe(unsigned int var) const
Accessor for side finite element object for scalar-valued variable var for the largest dimension in t...
 
void _do_elem_position_set(Real theta)
Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to th...
 
void side_values(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_values_vector) const
Fills a vector of values of the _system_vector at the all the quadrature points on the current elemen...
 
AlgebraicType _atype
Keep track of what type of algebra reinitialization is to be done.
 
unsigned char edge
Current edge for edge_* to examine.
 
std::map< FEType, std::unique_ptr< FEAbstract > > _edge_fe
 
Tensor side_hessian(unsigned int var, unsigned int qp) const
 
unsigned int get_mesh_x_var() const
Accessor for x-variable of moving mesh System.
 
void interior_values(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_values_vector) const
Fills a vector of values of the _system_vector at the all the quadrature points in the current elemen...
 
Tensor fixed_side_hessian(unsigned int var, unsigned int qp) const
 
void elem_position_set(Real theta)
Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to th...
 
void interior_hessians(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
Fills a vector of hessians of the _system_vector at the all the quadrature points in the current elem...
 
void interior_rate_gradient(unsigned int var, unsigned int qp, OutputType &u) const
 
void(FEMContext::* value_getter)(unsigned int, value_base *&, unsigned short) const
 
const NumericVector< Number > * _custom_solution
Data with which to do algebra reinitialization.
 
This class provides all data required for a physics package (e.g.
 
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
 
void attach_quadrature_rules()
Helper function for attaching quadrature rules.
 
void use_default_quadrature_rules(int extra_quadrature_order=0)
Use quadrature rules designed to over-integrate a mass matrix, plus extra_quadrature_order.
 
void _update_time_from_system(Real theta)
Update the time in the context object for the given value of theta, based on the values of "time" and...
 
virtual void set_mesh_system(System *sys)
Tells the FEMContext that system sys contains the isoparametric Lagrangian variables which correspond...
 
const typedef DenseSubVector< Number > &(DiffContext::* diff_subsolution_getter)(unsigned int) const
Helper typedef to simplify refactoring.
 
FEGenericBase< value_shape > value_base
 
Tensor point_hessian(unsigned int var, const Point &p) const
 
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
Reinitializes interior FE objects on the current geometric element.
 
void point_accel(unsigned int var, const Point &p, OutputType &u) const
 
virtual void edge_fe_reinit()
Reinitializes edge FE objects on the current geometric element.
 
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.
 
std::vector< std::unique_ptr< QBase > > _side_qrule
Quadrature rules for element sides The FEM context will try to find a quadrature rule that correctly ...
 
Gradient interior_gradient(unsigned int var, unsigned int qp) const
 
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
 
void set_mesh_y_var(unsigned int y_var)
Accessor for y-variable of moving mesh System.
 
void set_algebraic_type(const AlgebraicType atype)
Setting which determines whether to initialize algebraic structures (elem_*) on each element and set ...
 
TensorTools::MakeReal< OutputType >::type value_shape
 
FEGenericBase< hess_shape > hess_base
 
bool has_side_boundary_id(boundary_id_type id) const
Reports if the boundary id is found on the current side.
 
Number side_value(unsigned int var, unsigned int qp) const
 
AlgebraicType algebraic_type() const
 
const QBase & get_edge_qrule() const
Accessor for element edge quadrature rule.
 
This is the base class from which all geometric element types are derived.
 
const QBase & get_side_qrule() const
Accessor for element side quadrature rule for the dimension of the current _elem.
 
unsigned char get_dim() const
Accessor for cached mesh dimension.
 
AlgebraicType
Enum describing what data to use when initializing algebraic structures on each element.
 
Number fixed_interior_value(unsigned int var, unsigned int qp) const
 
FEType find_hardest_fe_type()
Helper function for creating quadrature rules.
 
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...
 
void some_value(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_value methods.
 
const std::set< unsigned char > & elem_dimensions() const
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
unsigned int get_mesh_z_var() const
Accessor for z-variable of moving mesh System.
 
bool _real_grad_fe_is_inf
 
virtual void elem_edge_reinit(Real theta) override
Resets the current time in the context.
 
unsigned char get_edge() const
Accessor for current edge of Elem object.
 
void interior_div(unsigned int var, unsigned int qp, OutputType &div_u) const
 
Tensor fixed_interior_hessian(unsigned int var, unsigned int qp) const
 
std::vector< std::vector< FEAbstract * > > _element_fe_var
Pointers to the same finite element objects, but indexed by variable number.
 
Number interior_value(unsigned int var, unsigned int qp) const
 
void set_mesh_z_var(unsigned int z_var)
Accessor for z-variable of moving mesh System.
 
Tensor fixed_point_hessian(unsigned int var, const Point &p) const
 
void ErrorVector unsigned int
 
unsigned char _dim
Cached dimension of largest dimension element in this mesh.
 
void set_mesh_x_var(unsigned int x_var)
Accessor for x-variable of moving mesh System.
 
std::vector< FEAbstract * > _edge_fe_var
 
virtual void elem_reinit(Real theta) override
Resets the current time in the context.
 
void side_hessians(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
Fills a vector of hessians of the _system_vector at the all the quadrature points on the current elem...
 
Gradient fixed_interior_gradient(unsigned int var, unsigned int qp) const
 
Elem & get_elem()
Accessor for current Elem object.
 
This class provides all data required for a physics package (e.g.