20 #ifndef LIBMESH_FEM_CONTEXT_H 21 #define LIBMESH_FEM_CONTEXT_H 30 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 43 template <
typename T>
class FEGenericBase;
44 typedef FEGenericBase<Real>
FEBase;
47 template <
typename T>
class NumericVector;
94 #ifdef LIBMESH_ENABLE_DEPRECATED 151 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 176 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 226 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 251 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 261 template<
typename OutputShape>
263 { this->get_element_fe<OutputShape>(var,fe,this->
get_dim()); }
280 template<
typename OutputShape>
282 unsigned short dim )
const;
298 template<
typename OutputShape>
300 { this->get_side_fe<OutputShape>(var,fe,this->
get_dim()); }
317 template<
typename OutputShape>
319 unsigned short dim )
const;
330 template<
typename OutputShape>
344 template<
typename OutputType>
347 OutputType & u)
const;
353 template<
typename OutputType>
356 std::vector<OutputType> & interior_values_vector)
const;
364 template<
typename OutputType>
367 OutputType & u)
const;
373 template<
typename OutputType>
376 std::vector<OutputType> & side_values_vector)
const;
387 template<
typename OutputType>
399 template<
typename OutputType>
401 OutputType & du)
const;
409 template<
typename OutputType>
412 std::vector<OutputType> & interior_gradients_vector)
const;
420 template<
typename OutputType>
423 OutputType & du)
const;
431 template<
typename OutputType>
434 std::vector<OutputType> & side_gradients_vector)
const;
445 template<
typename OutputType>
451 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 458 template<
typename OutputType>
461 OutputType & d2u)
const;
468 template<
typename OutputType>
471 std::vector<OutputType> & d2u_vals)
const;
479 template<
typename OutputType>
482 OutputType & d2u)
const;
489 template<
typename OutputType>
492 std::vector<OutputType> & d2u_vals)
const;
503 template<
typename OutputType>
509 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 516 template<
typename OutputType>
519 OutputType & u)
const;
525 template<
typename OutputType>
528 OutputType & u)
const;
534 template<
typename OutputType>
537 OutputType & u)
const;
544 template<
typename OutputType>
547 OutputType & u)
const;
553 template<
typename OutputType>
556 OutputType & u)
const;
562 template<
typename OutputType>
565 OutputType & u)
const;
573 template<
typename OutputType>
576 OutputType & u)
const;
584 template<
typename OutputType>
587 OutputType & u)
const;
598 template<
typename OutputType>
610 template<
typename OutputType>
613 OutputType & grad_u)
const;
621 template<
typename OutputType>
624 OutputType & grad_u)
const;
635 template<
typename OutputType>
641 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 648 template<
typename OutputType>
651 OutputType & hess_u)
const;
659 template<
typename OutputType>
662 OutputType & hess_u)
const;
673 template<
typename OutputType>
679 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 685 template<
typename OutputType>
688 OutputType & curl_u)
const;
697 template<
typename OutputType>
707 template<
typename OutputType>
710 OutputType & div_u)
const;
749 virtual void elem_fe_reinit(
const std::vector<Point> *
const pts =
nullptr);
866 {
return (this->
_elem !=
nullptr); }
872 { libmesh_assert(this->
_elem);
873 return *(this->
_elem); }
879 { libmesh_assert(this->
_elem);
880 return *(
const_cast<Elem *
>(this->
_elem)); }
900 {
return this->
_dim; }
989 template<
typename OutputShape>
1009 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 1014 template<
typename OutputShape>
1016 const FEType fe_type )
const;
1036 template <
typename OutputType>
1065 template<
typename OutputType,
1068 void some_value(
unsigned int var,
unsigned int qp, OutputType & u)
const;
1074 template<
typename OutputType,
1077 void some_gradient(
unsigned int var,
unsigned int qp, OutputType & u)
const;
1079 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 1084 template<
typename OutputType,
1087 void some_hessian(
unsigned int var,
unsigned int qp, OutputType & u)
const;
1095 std::vector<std::map<FEType, std::unique_ptr<FEAbstract>>>
_element_fe;
1096 std::vector<std::map<FEType, std::unique_ptr<FEAbstract>>>
_side_fe;
1097 std::map<FEType, std::unique_ptr<FEAbstract>>
_edge_fe;
1203 template<
typename OutputShape>
1206 unsigned short dim )
const 1210 fe = cast_ptr<FEGenericBase<OutputShape> *>( (
_element_fe_var[dim][var] ) );
1221 template<
typename OutputShape>
1224 unsigned short dim )
const 1227 libmesh_assert_less ( var, (
_side_fe_var[dim].size() ) );
1228 fe = cast_ptr<FEGenericBase<OutputShape> *>( (
_side_fe_var[dim][var] ) );
1235 libmesh_assert_less ( var, (
_side_fe_var[dim].size() ) );
1236 return cast_ptr<FEBase *>( (
_side_fe_var[dim][var] ) );
1239 template<
typename OutputShape>
1244 fe = cast_ptr<FEGenericBase<OutputShape> *>(
_edge_fe_var[var] );
1257 #endif // LIBMESH_FEM_CONTEXT_H void point_accel(unsigned int var, const Point &p, OutputType &u) const
Manages the family, order, etc. parameters for a given FE.
unsigned char get_edge() const
virtual void nonlocal_reinit(Real theta) override
void(FEMContext::* value_getter)(unsigned int, value_base *&, unsigned short) const
Number fixed_side_value(unsigned int var, unsigned int qp) const
Number side_value(unsigned int var, unsigned int qp) const
Tensor fixed_point_hessian(unsigned int var, const Point &p) const
void set_mesh_z_var(unsigned int z_var)
void point_rate(unsigned int var, const Point &p, OutputType &u) const
Number interior_value(unsigned int var, unsigned int qp) const
void interior_div(unsigned int var, unsigned int qp, OutputType &div_u) const
void elem_position_set(Real theta)
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
virtual void pre_fe_reinit(const System &, const Elem *e)
unsigned int get_mesh_x_var() const
void set_mesh_x_var(unsigned int x_var)
const NumericVector< Number > * _custom_solution
void set_elem(const Elem *e)
void get_edge_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
std::unique_ptr< FEGenericBase< RealGradient > > _real_grad_fe
unsigned char get_elem_dim() const
const Elem & get_elem() const
bool _real_grad_fe_is_inf
TensorTools::MakeReal< OutputType >::type value_shape
const QBase & get_element_qrule(unsigned short dim) const
const DenseSubVector< Number > &(DiffContext::* diff_subsolution_getter)(unsigned int) const
FEGenericBase< value_shape > value_base
virtual void elem_edge_reinit(Real theta) override
The base class for all geometric element types.
void set_algebraic_type(const AlgebraicType atype)
int _extra_quadrature_order
static const Real TOLERANCE
FEGenericBase< hess_shape > hess_base
void interior_gradients(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_gradients_vector) const
bool has_side_boundary_id(boundary_id_type id) const
unsigned char get_side() const
Number point_value(unsigned int var, const Point &p) const
DiffContext(const System &)
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _element_fe
unsigned char get_dim() const
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
std::unique_ptr< FEGenericBase< Real > > _real_fe
Number fixed_interior_value(unsigned int var, unsigned int qp) const
virtual void set_mesh_system(System *sys)
Tensor fixed_side_hessian(unsigned int var, unsigned int qp) const
std::set< unsigned char > _elem_dims
Gradient interior_gradient(unsigned int var, unsigned int qp) const
void point_curl(unsigned int var, const Point &p, OutputType &curl_u, const Real tolerance=TOLERANCE) const
void side_values(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_values_vector) const
Tensor side_hessian(unsigned int var, unsigned int qp) const
virtual void side_fe_reinit()
void interior_hessians(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
std::map< FEType, std::unique_ptr< FEAbstract > > _edge_fe
void _update_time_from_system(Real theta)
void init_internal_data(const System &sys)
Gradient fixed_point_gradient(unsigned int var, const Point &p) const
const QBase & get_side_qrule(unsigned short dim) const
void interior_values(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_values_vector) const
Manages consistently variables, degrees of freedom, and coefficient vectors.
void some_value(unsigned int var, unsigned int qp, OutputType &u) const
void some_hessian(unsigned int var, unsigned int qp, OutputType &u) const
Tensor interior_hessian(unsigned int var, unsigned int qp) const
Tensor point_hessian(unsigned int var, const Point &p) const
Gradient fixed_interior_gradient(unsigned int var, unsigned int qp) const
Used by the Mesh to keep track of boundary nodes and elements.
TensorTools::MakeReal< Rank1Decrement >::type grad_shape
TensorTools::MakeReal< Rank2Decrement >::type hess_shape
const BoundaryInfo & _boundary_info
std::vector< std::vector< FEAbstract * > > _element_fe_var
FEGenericBase< grad_shape > grad_base
void(FEMContext::* hess_getter)(unsigned int, hess_base *&, unsigned short) const
FEBase * get_side_fe(unsigned int var) const
std::vector< boundary_id_type > side_boundary_ids() const
FEGenericBase< Real > FEBase
TensorTools::DecrementRank< Rank1Decrement >::type Rank2Decrement
virtual void elem_side_reinit(Real theta) override
Number fixed_point_value(unsigned int var, const Point &p) const
virtual void edge_fe_reinit()
std::vector< FEAbstract * > _edge_fe_var
void interior_rate(unsigned int var, unsigned int qp, OutputType &u) const
std::vector< std::unique_ptr< QBase > > _side_qrule
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _side_fe
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< std::unique_ptr< QBase > > _element_qrule
void side_accel(unsigned int var, unsigned int qp, OutputType &u) const
void interior_curl(unsigned int var, unsigned int qp, OutputType &curl_u) const
TensorTools::DecrementRank< OutputType >::type Rank1Decrement
void set_mesh_y_var(unsigned int y_var)
std::unique_ptr< QBase > _edge_qrule
FEMContext(const System &sys)
const System * get_mesh_system() const
void side_gradients(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_gradients_vector) const
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Gradient fixed_side_gradient(unsigned int var, unsigned int qp) const
unsigned int get_mesh_y_var() const
void(FEMContext::* grad_getter)(unsigned int, grad_base *&, unsigned short) const
AlgebraicType algebraic_type() const
const QBase & get_edge_qrule() const
void some_gradient(unsigned int var, unsigned int qp, OutputType &u) const
void set_custom_solution(const NumericVector< Number > *custom_sol)
std::vector< std::vector< FEAbstract * > > _side_fe_var
FEGenericBase< OutputShape > * build_new_fe(const FEGenericBase< OutputShape > *fe, const Point &p, const Real tolerance=TOLERANCE) const
Gradient point_gradient(unsigned int var, const Point &p) const
FEGenericBase< OutputShape > * cached_fe(const unsigned int elem_dim, const FEType fe_type) const
void _do_elem_position_set(Real theta)
Gradient side_gradient(unsigned int var, unsigned int qp) const
FEBase * get_element_fe(unsigned int var) const
void side_hessians(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
void side_rate(unsigned int var, unsigned int qp, OutputType &u) const
A geometric point in (x,y,z) space.
System * get_mesh_system()
Base class for all quadrature families and orders.
const QBase & get_element_qrule() const
virtual void elem_reinit(Real theta) override
unsigned int get_mesh_z_var() const
const std::set< unsigned char > & elem_dimensions() const
void interior_accel(unsigned int var, unsigned int qp, OutputType &u) const
const QBase & get_side_qrule() const
Tensor fixed_interior_hessian(unsigned int var, unsigned int qp) const