libMesh
fem_context.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_FEM_CONTEXT_H
21 #define LIBMESH_FEM_CONTEXT_H
22 
23 // Local Includes
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"
29 
30 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
31 #include "libmesh/tensor_value.h"
32 #endif
33 
34 // C++ includes
35 #include <map>
36 #include <set>
37 
38 namespace libMesh
39 {
40 
41 // Forward Declarations
42 class BoundaryInfo;
43 class Elem;
44 template <typename T> class FEGenericBase;
45 typedef FEGenericBase<Real> FEBase;
46 class QBase;
47 class Point;
48 template <typename T> class NumericVector;
49 
62 class FEMContext : public DiffContext
63 {
64 public:
65 
69  explicit
70  FEMContext (const System & sys);
71 
76  explicit
77  FEMContext (const System & sys, int extra_quadrature_order);
78 
82  virtual ~FEMContext ();
83 
88  void use_default_quadrature_rules(int extra_quadrature_order=0);
89 
94  void use_unweighted_quadrature_rules(int extra_quadrature_order=0);
95 
100 
107 #ifdef LIBMESH_ENABLE_DEPRECATED
108  std::vector<boundary_id_type> side_boundary_ids() const;
109 #endif
110 
114  void side_boundary_ids(std::vector<boundary_id_type> & vec_to_fill) const;
115 
122  Number interior_value(unsigned int var, unsigned int qp) const;
123 
130  Number side_value(unsigned int var, unsigned int qp) const;
131 
138  Number point_value(unsigned int var, const Point & p) const;
139 
146  Gradient interior_gradient(unsigned int var, unsigned int qp) const;
147 
154  Gradient side_gradient(unsigned int var, unsigned int qp) const;
155 
162  Gradient point_gradient(unsigned int var, const Point & p) const;
163 
164 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
165 
171  Tensor interior_hessian(unsigned int var, unsigned int qp) const;
172 
179  Tensor side_hessian(unsigned int var, unsigned int qp) const;
180 
187  Tensor point_hessian(unsigned int var, const Point & p) const;
188 
189 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
190 
197  Number fixed_interior_value(unsigned int var, unsigned int qp) const;
198 
205  Number fixed_side_value(unsigned int var, unsigned int qp) const;
206 
213  Number fixed_point_value(unsigned int var, const Point & p) const;
214 
221  Gradient fixed_interior_gradient(unsigned int var, unsigned int qp) const;
222 
229  Gradient fixed_side_gradient(unsigned int var, unsigned int qp) const;
230 
237  Gradient fixed_point_gradient(unsigned int var, const Point & p) const;
238 
239 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
240 
246  Tensor fixed_interior_hessian(unsigned int var, unsigned int qp) const;
247 
254  Tensor fixed_side_hessian(unsigned int var, unsigned int qp) const;
255 
262  Tensor fixed_point_hessian (unsigned int var, const Point & p) const;
263 
264 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
265 
274  template<typename OutputShape>
275  void get_element_fe( unsigned int var, FEGenericBase<OutputShape> *& fe ) const
276  { this->get_element_fe<OutputShape>(var,fe,this->get_dim()); }
277 
286  FEBase * get_element_fe( unsigned int var ) const
287  { return this->get_element_fe(var,this->get_dim()); }
288 
293  template<typename OutputShape>
294  void get_element_fe( unsigned int var, FEGenericBase<OutputShape> *& fe,
295  unsigned short dim ) const;
296 
301  FEBase * get_element_fe( unsigned int var, unsigned short dim ) const;
302 
311  template<typename OutputShape>
312  void get_side_fe( unsigned int var, FEGenericBase<OutputShape> *& fe ) const
313  { this->get_side_fe<OutputShape>(var,fe,this->get_dim()); }
314 
323  FEBase * get_side_fe( unsigned int var ) const
324  { return this->get_side_fe(var,this->get_dim()); }
325 
330  template<typename OutputShape>
331  void get_side_fe( unsigned int var, FEGenericBase<OutputShape> *& fe,
332  unsigned short dim ) const;
333 
338  FEBase * get_side_fe( unsigned int var, unsigned short dim ) const;
339 
343  template<typename OutputShape>
344  void get_edge_fe( unsigned int var, FEGenericBase<OutputShape> *& fe ) const;
345 
349  FEBase * get_edge_fe( unsigned int var ) const;
350 
357  template<typename OutputType>
358  void interior_value(unsigned int var,
359  unsigned int qp,
360  OutputType & u) const;
361 
366  template<typename OutputType>
367  void interior_values(unsigned int var,
368  const NumericVector<Number> & _system_vector,
369  std::vector<OutputType> & interior_values_vector) const;
370 
377  template<typename OutputType>
378  void side_value(unsigned int var,
379  unsigned int qp,
380  OutputType & u) const;
381 
386  template<typename OutputType>
387  void side_values(unsigned int var,
388  const NumericVector<Number> & _system_vector,
389  std::vector<OutputType> & side_values_vector) const;
390 
400  template<typename OutputType>
401  void point_value(unsigned int var,
402  const Point & p,
403  OutputType & u,
404  const Real tolerance = TOLERANCE) const;
405 
412  template<typename OutputType>
413  void interior_gradient(unsigned int var, unsigned int qp,
414  OutputType & du) const;
415 
422  template<typename OutputType>
423  void interior_gradients(unsigned int var,
424  const NumericVector<Number> & _system_vector,
425  std::vector<OutputType> & interior_gradients_vector) const;
426 
433  template<typename OutputType>
434  void side_gradient(unsigned int var,
435  unsigned int qp,
436  OutputType & du) const;
437 
444  template<typename OutputType>
445  void side_gradients(unsigned int var,
446  const NumericVector<Number> & _system_vector,
447  std::vector<OutputType> & side_gradients_vector) const;
448 
458  template<typename OutputType>
459  void point_gradient(unsigned int var,
460  const Point & p,
461  OutputType & grad_u,
462  const Real tolerance = TOLERANCE) const;
463 
464 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
465 
471  template<typename OutputType>
472  void interior_hessian(unsigned int var,
473  unsigned int qp,
474  OutputType & d2u) const;
475 
481  template<typename OutputType>
482  void interior_hessians(unsigned int var,
483  const NumericVector<Number> & _system_vector,
484  std::vector<OutputType> & d2u_vals) const;
485 
492  template<typename OutputType>
493  void side_hessian(unsigned int var,
494  unsigned int qp,
495  OutputType & d2u) const;
496 
502  template<typename OutputType>
503  void side_hessians(unsigned int var,
504  const NumericVector<Number> & _system_vector,
505  std::vector<OutputType> & d2u_vals) const;
506 
516  template<typename OutputType>
517  void point_hessian(unsigned int var,
518  const Point & p,
519  OutputType & hess_u,
520  const Real tolerance = TOLERANCE) const;
521 
522 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
523 
529  template<typename OutputType>
530  void interior_rate(unsigned int var,
531  unsigned int qp,
532  OutputType & u) const;
533 
534 
540  template<typename OutputType>
541  void interior_rate_gradient(unsigned int var,
542  unsigned int qp,
543  OutputType & u) const;
544 
545 
550  template<typename OutputType>
551  void side_rate(unsigned int var,
552  unsigned int qp,
553  OutputType & u) const;
554 
559  template<typename OutputType>
560  void point_rate(unsigned int var,
561  const Point & p,
562  OutputType & u) const;
563 
569  template<typename OutputType>
570  void interior_accel(unsigned int var,
571  unsigned int qp,
572  OutputType & u) const;
573 
578  template<typename OutputType>
579  void side_accel(unsigned int var,
580  unsigned int qp,
581  OutputType & u) const;
582 
587  template<typename OutputType>
588  void point_accel(unsigned int var,
589  const Point & p,
590  OutputType & u) const;
591 
598  template<typename OutputType>
599  void fixed_interior_value(unsigned int var,
600  unsigned int qp,
601  OutputType & u) const;
602 
609  template<typename OutputType>
610  void fixed_side_value(unsigned int var,
611  unsigned int qp,
612  OutputType & u) const;
613 
623  template<typename OutputType>
624  void fixed_point_value(unsigned int var,
625  const Point & p,
626  OutputType & u,
627  const Real tolerance = TOLERANCE) const;
628 
635  template<typename OutputType>
636  void fixed_interior_gradient(unsigned int var,
637  unsigned int qp,
638  OutputType & grad_u) const;
639 
646  template<typename OutputType>
647  void fixed_side_gradient(unsigned int var,
648  unsigned int qp,
649  OutputType & grad_u) const;
650 
660  template<typename OutputType>
661  void fixed_point_gradient(unsigned int var,
662  const Point & p,
663  OutputType & grad_u,
664  const Real tolerance = TOLERANCE) const;
665 
666 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
667 
673  template<typename OutputType>
674  void fixed_interior_hessian(unsigned int var,
675  unsigned int qp,
676  OutputType & hess_u) const;
677 
684  template<typename OutputType>
685  void fixed_side_hessian(unsigned int var,
686  unsigned int qp,
687  OutputType & hess_u) const;
688 
698  template<typename OutputType>
699  void fixed_point_hessian(unsigned int var,
700  const Point & p,
701  OutputType & hess_u,
702  const Real tolerance = TOLERANCE) const;
703 
704 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
705 
710  template<typename OutputType>
711  void interior_curl(unsigned int var,
712  unsigned int qp,
713  OutputType & curl_u) const;
714 
722  template<typename OutputType>
723  void point_curl(unsigned int var,
724  const Point & p,
725  OutputType & curl_u,
726  const Real tolerance = TOLERANCE) const;
727 
732  template<typename OutputType>
733  void interior_div(unsigned int var,
734  unsigned int qp,
735  OutputType & div_u) const;
736 
737  // should be protected:
743  virtual void elem_reinit(Real theta) override;
744 
750  virtual void elem_side_reinit(Real theta) override;
751 
757  virtual void elem_edge_reinit(Real theta) override;
758 
763  virtual void nonlocal_reinit(Real theta) override;
764 
768  virtual void pre_fe_reinit(const System &,
769  const Elem * e);
770 
774  virtual void elem_fe_reinit(const std::vector<Point> * const pts = nullptr);
775 
779  virtual void side_fe_reinit();
780 
784  virtual void edge_fe_reinit();
785 
790  const QBase & get_element_qrule() const
791  { return this->get_element_qrule(this->get_elem_dim()); }
792 
797  const QBase & get_side_qrule() const
798  { return this->get_side_qrule(this->get_elem_dim()); }
799 
803  const QBase & get_element_qrule( unsigned short dim ) const
805  return *(this->_element_qrule[dim]); }
806 
810  const QBase & get_side_qrule( unsigned short dim ) const
811  {
813  return *(this->_side_qrule[dim]);
814  }
815 
819  const QBase & get_edge_qrule() const
820  { return *(this->_edge_qrule); }
821 
830  virtual void set_mesh_system(System * sys)
831  { this->_mesh_sys = sys; }
832 
836  const System * get_mesh_system() const
837  { return this->_mesh_sys; }
838 
843  { return this->_mesh_sys; }
844 
848  unsigned int get_mesh_x_var() const
849  { return _mesh_x_var; }
850 
856  void set_mesh_x_var(unsigned int x_var)
857  { _mesh_x_var = x_var; }
858 
862  unsigned int get_mesh_y_var() const
863  { return _mesh_y_var; }
864 
870  void set_mesh_y_var(unsigned int y_var)
871  { _mesh_y_var = y_var; }
872 
876  unsigned int get_mesh_z_var() const
877  { return _mesh_z_var; }
878 
884  void set_mesh_z_var(unsigned int z_var)
885  { _mesh_z_var = z_var; }
886 
890  bool has_elem() const
891  { return (this->_elem != nullptr); }
892 
896  const Elem & get_elem() const
897  { libmesh_assert(this->_elem);
898  return *(this->_elem); }
899 
904  { libmesh_assert(this->_elem);
905  return *(const_cast<Elem *>(this->_elem)); }
906 
910  unsigned char get_side() const
911  { return side; }
912 
916  unsigned char get_edge() const
917  { return edge; }
918 
924  unsigned char get_dim() const
925  { return this->_dim; }
926 
931  unsigned char get_elem_dim() const
932  { return _elem_dim; }
933 
938  const std::set<unsigned char> & elem_dimensions() const
939  { return _elem_dims; }
940 
946  void elem_position_set(Real theta);
947 
952  void elem_position_get();
953 
958  enum AlgebraicType { NONE = 0, // Do not reinitialize dof_indices
959  DOFS_ONLY, // Reinitialize dof_indices, not
960  // algebraic structures
961  CURRENT, // Use dof_indices, current solution
962  OLD, // Use old_dof_indices, custom solution
963  OLD_DOFS_ONLY}; // Reinitialize old_dof_indices, not
964  // algebraic structures
965 
974  { _atype = atype; }
975 
976  /*
977  * Get the current AlgebraicType setting
978  */
979  AlgebraicType algebraic_type() const { return _atype; }
980 
989  { _custom_solution = custom_sol; }
990 
995  void set_jacobian_tolerance(Real tol);
996 
1001 
1006 
1010  unsigned char side;
1011 
1015  unsigned char edge;
1016 
1021 
1025  void attach_quadrature_rules();
1026 
1030  template<typename OutputShape>
1032  const Point & p,
1033  const Real tolerance = TOLERANCE) const;
1034 
1035 protected:
1036 
1041 
1046 
1047  mutable std::unique_ptr<FEGenericBase<Real>> _real_fe;
1048  mutable std::unique_ptr<FEGenericBase<RealGradient>> _real_grad_fe;
1049 
1050 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1051  mutable bool _real_fe_is_inf;
1052  mutable bool _real_grad_fe_is_inf;
1053 #endif
1054 
1055  template<typename OutputShape>
1056  FEGenericBase<OutputShape> * cached_fe( const unsigned int elem_dim,
1057  const FEType fe_type ) const;
1058 
1062  void set_elem( const Elem * e );
1063 
1064  // gcc-3.4, oracle 12.3 require this typedef to be public
1065  // in order to use it in a return type
1066 public:
1067 
1071  typedef const DenseSubVector<Number> & (DiffContext::*diff_subsolution_getter)(unsigned int) const;
1072 
1073 protected:
1077  template <typename OutputType>
1078  struct FENeeded
1079  {
1080  // Rank decrementer helper types
1083 
1084  // Typedefs for "Value getter" function pointer
1087  typedef void (FEMContext::*value_getter) (unsigned int, value_base *&, unsigned short) const;
1088 
1089  // Typedefs for "Grad getter" function pointer
1092  typedef void (FEMContext::*grad_getter) (unsigned int, grad_base *&, unsigned short) const;
1093 
1094  // Typedefs for "Hessian getter" function pointer
1097  typedef void (FEMContext::*hess_getter) (unsigned int, hess_base *&, unsigned short) const;
1098  };
1099 
1100 
1101 
1106  template<typename OutputType,
1107  typename FENeeded<OutputType>::value_getter fe_getter,
1108  diff_subsolution_getter subsolution_getter>
1109  void some_value(unsigned int var, unsigned int qp, OutputType & u) const;
1110 
1115  template<typename OutputType,
1116  typename FENeeded<OutputType>::grad_getter fe_getter,
1117  diff_subsolution_getter subsolution_getter>
1118  void some_gradient(unsigned int var, unsigned int qp, OutputType & u) const;
1119 
1120 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1121 
1125  template<typename OutputType,
1126  typename FENeeded<OutputType>::hess_getter fe_getter,
1127  diff_subsolution_getter subsolution_getter>
1128  void some_hessian(unsigned int var, unsigned int qp, OutputType & u) const;
1129 #endif
1130 
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;
1139 
1140 
1147  std::vector<std::vector<FEAbstract *>> _element_fe_var;
1148  std::vector<std::vector<FEAbstract *>> _side_fe_var;
1149  std::vector<FEAbstract *> _edge_fe_var;
1150 
1156 
1160  const Elem * _elem;
1161 
1165  unsigned char _dim;
1166 
1170  unsigned char _elem_dim;
1171 
1176  std::set<unsigned char> _elem_dims;
1177 
1184  std::vector<std::unique_ptr<QBase>> _element_qrule;
1185 
1192  std::vector<std::unique_ptr<QBase>> _side_qrule;
1193 
1201  std::unique_ptr<QBase> _edge_qrule;
1202 
1207 
1208 private:
1212  void init_internal_data(const System & sys);
1213 
1222  void _do_elem_position_set(Real theta);
1223 
1229  void _update_time_from_system(Real theta);
1230 };
1231 
1232 
1233 
1234 // ------------------------------------------------------------
1235 // FEMContext inline methods
1236 
1237 inline
1239 {
1240  if (_mesh_sys)
1241  this->_do_elem_position_set(theta);
1242 }
1243 
1244 template<typename OutputShape>
1245 inline
1247  unsigned short dim ) const
1248 {
1249  libmesh_assert( !_element_fe_var[dim].empty() );
1250  libmesh_assert_less ( var, (_element_fe_var[dim].size() ) );
1251  fe = cast_ptr<FEGenericBase<OutputShape> *>( (_element_fe_var[dim][var] ) );
1252 }
1253 
1254 inline
1255 FEBase * FEMContext::get_element_fe( unsigned int var, unsigned short dim ) const
1256 {
1257  libmesh_assert( !_element_fe_var[dim].empty() );
1258  libmesh_assert_less ( var, (_element_fe_var[dim].size() ) );
1259  return cast_ptr<FEBase *>( (_element_fe_var[dim][var] ) );
1260 }
1261 
1262 template<typename OutputShape>
1263 inline
1265  unsigned short dim ) const
1266 {
1267  libmesh_assert( !_side_fe_var[dim].empty() );
1268  libmesh_assert_less ( var, (_side_fe_var[dim].size() ) );
1269  fe = cast_ptr<FEGenericBase<OutputShape> *>( (_side_fe_var[dim][var] ) );
1270 }
1271 
1272 inline
1273 FEBase * FEMContext::get_side_fe( unsigned int var, unsigned short dim ) const
1274 {
1275  libmesh_assert( !_side_fe_var[dim].empty() );
1276  libmesh_assert_less ( var, (_side_fe_var[dim].size() ) );
1277  return cast_ptr<FEBase *>( (_side_fe_var[dim][var] ) );
1278 }
1279 
1280 template<typename OutputShape>
1281 inline
1282 void FEMContext::get_edge_fe( unsigned int var, FEGenericBase<OutputShape> *& fe ) const
1283 {
1284  libmesh_assert_less ( var, _edge_fe_var.size() );
1285  fe = cast_ptr<FEGenericBase<OutputShape> *>( _edge_fe_var[var] );
1286 }
1287 
1288 inline
1289 FEBase * FEMContext::get_edge_fe( unsigned int var ) const
1290 {
1291  libmesh_assert_less ( var, _edge_fe_var.size() );
1292  return cast_ptr<FEBase *>( _edge_fe_var[var] );
1293 }
1294 
1295 
1296 } // namespace libMesh
1297 
1298 #endif // LIBMESH_FEM_CONTEXT_H
libMesh::QBase
The QBase class provides the basic functionality from which various quadrature rules can be derived.
Definition: quadrature.h:61
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::FEMContext::FENeeded::hess_shape
TensorTools::MakeReal< Rank2Decrement >::type hess_shape
Definition: fem_context.h:1095
libMesh::FEMContext::FENeeded
Helper nested class for C++03-compatible "template typedef".
Definition: fem_context.h:1078
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::FEMContext::init_internal_data
void init_internal_data(const System &sys)
Helper function used in constructors to set up internal data.
Definition: fem_context.C:176
libMesh::FEMContext::get_element_qrule
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:790
libMesh::FEMContext::fixed_point_gradient
Gradient fixed_point_gradient(unsigned int var, const Point &p) const
Definition: fem_context.C:1214
libMesh::DiffContext::DiffContext
DiffContext(const System &)
Constructor.
Definition: diff_context.C:30
libMesh::BoundaryInfo
The BoundaryInfo class contains information relevant to boundary conditions including storing faces,...
Definition: boundary_info.h:57
libMesh::FEMContext::side_gradient
Gradient side_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:678
libMesh::FEMContext::side
unsigned char side
Current side for side_* to examine.
Definition: fem_context.h:1010
libMesh::FEMContext::elem_position_get
void elem_position_get()
Uses the geometry of elem to set the coordinate data specified by mesh_*_position configuration.
Definition: fem_context.C:1485
libMesh::FEMContext::FENeeded::grad_base
FEGenericBase< grad_shape > grad_base
Definition: fem_context.h:1091
libMesh::FEMContext::get_element_qrule
const QBase & get_element_qrule(unsigned short dim) const
Accessor for element interior quadrature rule.
Definition: fem_context.h:803
libMesh::FEMContext::_real_grad_fe
std::unique_ptr< FEGenericBase< RealGradient > > _real_grad_fe
Definition: fem_context.h:1048
libMesh::FEMContext::FENeeded::hess_getter
void(FEMContext::* hess_getter)(unsigned int, hess_base *&, unsigned short) const
Definition: fem_context.h:1097
libMesh::FEMContext::~FEMContext
virtual ~FEMContext()
Destructor.
Definition: fem_context.C:248
libMesh::FEMContext::nonlocal_reinit
virtual void nonlocal_reinit(Real theta) override
Gives derived classes the opportunity to reinitialize data needed for nonlocal calculations at a new ...
Definition: fem_context.C:1426
libMesh::FEMContext::FEMContext
FEMContext(const System &sys)
Constructor.
Definition: fem_context.C:37
libMesh::FEMContext::_element_qrule
std::vector< std::unique_ptr< QBase > > _element_qrule
Quadrature rule for element interior.
Definition: fem_context.h:1184
libMesh::FEMContext::FENeeded::grad_shape
TensorTools::MakeReal< Rank1Decrement >::type grad_shape
Definition: fem_context.h:1090
libMesh::FEMContext::side_accel
void side_accel(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1361
libMesh::FEMContext::get_mesh_y_var
unsigned int get_mesh_y_var() const
Accessor for y-variable of moving mesh System.
Definition: fem_context.h:862
libMesh::FEMContext::elem_side_reinit
virtual void elem_side_reinit(Real theta) override
Resets the current time in the context.
Definition: fem_context.C:1396
libMesh::FEMContext::_edge_qrule
std::unique_ptr< QBase > _edge_qrule
Quadrature rules for element edges.
Definition: fem_context.h:1201
libMesh::FEMContext::side_boundary_ids
std::vector< boundary_id_type > side_boundary_ids() const
Lists the boundary ids found on the current side.
Definition: fem_context.C:261
libMesh::FEMContext::side_gradients
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...
Definition: fem_context.C:723
libMesh::FEMContext::fixed_side_gradient
Gradient fixed_side_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1118
libMesh::FEMContext::side_fe_reinit
virtual void side_fe_reinit()
Reinitializes side FE objects on the current geometric element.
Definition: fem_context.C:1457
libMesh::FEMContext::FENeeded::Rank2Decrement
TensorTools::DecrementRank< Rank1Decrement >::type Rank2Decrement
Definition: fem_context.h:1082
libMesh::FEMContext::build_new_fe
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.
Definition: fem_context.C:1952
libMesh::FEMContext::set_elem
void set_elem(const Elem *e)
Helper function to promote accessor usage.
Definition: fem_context.C:1855
libMesh::FEMContext::interior_gradients
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...
Definition: fem_context.C:453
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::FEMContext::get_side
unsigned char get_side() const
Accessor for current side of Elem object.
Definition: fem_context.h:910
libMesh::FEMContext::pre_fe_reinit
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
Definition: fem_context.C:1642
libMesh::FEGenericBase
This class forms the foundation from which generic finite elements may be derived.
Definition: exact_error_estimator.h:39
libMesh::FEMContext::_side_fe
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _side_fe
Definition: fem_context.h:1137
libMesh::FEMContext::get_element_fe
FEBase * get_element_fe(unsigned int var) const
Accessor for interior finite element object for scalar-valued variable var for the largest dimension ...
Definition: fem_context.h:286
libMesh::FEMContext::get_elem_dim
unsigned char get_elem_dim() const
Definition: fem_context.h:931
libMesh::FEMContext::get_elem
const Elem & get_elem() const
Accessor for current Elem object.
Definition: fem_context.h:896
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::FEMContext::get_edge_fe
void get_edge_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge (3D only!) finite element object for variable var.
Definition: fem_context.h:1282
libMesh::FEMContext::OLD_DOFS_ONLY
Definition: fem_context.h:963
libMesh::FEMContext::get_mesh_system
System * get_mesh_system()
Accessor for moving mesh System.
Definition: fem_context.h:842
libMesh::FEMContext::interior_accel
void interior_accel(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1349
libMesh::FEMContext::interior_curl
void interior_curl(unsigned int var, unsigned int qp, OutputType &curl_u) const
Definition: fem_context.C:557
libMesh::FEMContext::FENeeded::grad_getter
void(FEMContext::* grad_getter)(unsigned int, grad_base *&, unsigned short) const
Definition: fem_context.h:1092
libMesh::FEMContext::FENeeded::Rank1Decrement
TensorTools::DecrementRank< OutputType >::type Rank1Decrement
Definition: fem_context.h:1081
libMesh::boundary_id_type
int8_t boundary_id_type
Definition: id_types.h:51
libMesh::FEMContext::_boundary_info
const BoundaryInfo & _boundary_info
Saved reference to BoundaryInfo on the mesh for this System.
Definition: fem_context.h:1155
libMesh::TensorTools::MakeReal::type
T type
Definition: tensor_tools.h:250
libMesh::TensorTools::DecrementRank::type
T type
Definition: tensor_tools.h:158
libMesh::FEMContext::point_value
Number point_value(unsigned int var, const Point &p) const
Definition: fem_context.C:833
libMesh::FEMContext::set_custom_solution
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.
Definition: fem_context.h:988
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::FEMContext::_extra_quadrature_order
int _extra_quadrature_order
The extra quadrature order for this context.
Definition: fem_context.h:1206
libMesh::FEMContext::fixed_point_value
Number fixed_point_value(unsigned int var, const Point &p) const
Definition: fem_context.C:1168
libMesh::FEMContext::_elem_dim
unsigned char _elem_dim
Cached dimension of this->_elem.
Definition: fem_context.h:1170
libMesh::VectorValue< Number >
libMesh::FEMContext::some_gradient
void some_gradient(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_gradient methods.
Definition: fem_context.C:308
libMesh::NumericVector< Number >
libMesh::FEMContext::_side_fe_var
std::vector< std::vector< FEAbstract * > > _side_fe_var
Definition: fem_context.h:1148
libMesh::FEMContext::use_unweighted_quadrature_rules
void use_unweighted_quadrature_rules(int extra_quadrature_order=0)
Use quadrature rules designed to exactly integrate unweighted undistorted basis functions,...
Definition: fem_context.C:155
libMesh::FEMContext::fixed_side_value
Number fixed_side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1094
libMesh::FEMContext::point_rate
void point_rate(unsigned int var, const Point &p, OutputType &u) const
libMesh::FEMContext::cached_fe
FEGenericBase< OutputShape > * cached_fe(const unsigned int elem_dim, const FEType fe_type) const
libMesh::FEMContext::point_curl
void point_curl(unsigned int var, const Point &p, OutputType &curl_u, const Real tolerance=TOLERANCE) const
Definition: fem_context.C:982
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::FEMContext::interior_rate
void interior_rate(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1318
libMesh::FEMContext::_real_fe_is_inf
bool _real_fe_is_inf
Definition: fem_context.h:1051
libMesh::FEMContext::_elem_dims
std::set< unsigned char > _elem_dims
Cached dimensions of elements in the mesh, plus dimension 0 if SCALAR variables are in use.
Definition: fem_context.h:1176
libMesh::FEMContext::some_hessian
void some_hessian(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_hessian methods.
Definition: fem_context.C:341
libMesh::FEMContext::side_rate
void side_rate(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1339
libMesh::FEMContext::_elem
const Elem * _elem
Current element for element_* to examine.
Definition: fem_context.h:1160
libMesh::FEMContext::set_jacobian_tolerance
void set_jacobian_tolerance(Real tol)
Calls set_jacobian_tolerance() on all the FE objects controlled by this class.
Definition: fem_context.C:1545
libMesh::FEMContext::_mesh_sys
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: fem_context.h:1000
libMesh::FEMContext::has_elem
bool has_elem() const
Test for current Elem object.
Definition: fem_context.h:890
libMesh::FEBase
FEGenericBase< Real > FEBase
Definition: exact_error_estimator.h:39
libMesh::FEMContext::get_side_qrule
const QBase & get_side_qrule(unsigned short dim) const
Accessor for element side quadrature rule.
Definition: fem_context.h:810
libMesh::FEMContext::interior_hessian
Tensor interior_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:491
libMesh::FEMContext::_real_fe
std::unique_ptr< FEGenericBase< Real > > _real_fe
Definition: fem_context.h:1047
libMesh::FEMContext::get_mesh_system
const System * get_mesh_system() const
Accessor for moving mesh System.
Definition: fem_context.h:836
libMesh::TensorValue
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:53
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::FEMContext::_element_fe
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _element_fe
Finite element objects for each variable's interior, sides and edges.
Definition: fem_context.h:1136
libMesh::FEMContext::_mesh_y_var
unsigned int _mesh_y_var
Definition: fem_context.h:1005
libMesh::FEMContext::point_gradient
Gradient point_gradient(unsigned int var, const Point &p) const
Definition: fem_context.C:879
libMesh::FEMContext::get_side_fe
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...
Definition: fem_context.h:323
libMesh::FEMContext::_do_elem_position_set
void _do_elem_position_set(Real theta)
Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to th...
Definition: fem_context.C:1566
libMesh::FEMContext::side_values
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...
Definition: fem_context.C:643
libMesh::FEMContext::_atype
AlgebraicType _atype
Keep track of what type of algebra reinitialization is to be done.
Definition: fem_context.h:1040
libMesh::FEMContext::edge
unsigned char edge
Current edge for edge_* to examine.
Definition: fem_context.h:1015
libMesh::FEMContext::_edge_fe
std::map< FEType, std::unique_ptr< FEAbstract > > _edge_fe
Definition: fem_context.h:1138
libMesh::FEMContext::side_hessian
Tensor side_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:761
libMesh::FEMContext::get_mesh_x_var
unsigned int get_mesh_x_var() const
Accessor for x-variable of moving mesh System.
Definition: fem_context.h:848
libMesh::FEMContext::interior_values
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...
Definition: fem_context.C:391
libMesh::FEMContext::fixed_side_hessian
Tensor fixed_side_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1143
libMesh::FEMContext::elem_position_set
void elem_position_set(Real theta)
Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to th...
Definition: fem_context.h:1238
libMesh::FEMContext::interior_hessians
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...
Definition: fem_context.C:515
libMesh::FEMContext::interior_rate_gradient
void interior_rate_gradient(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1328
libMesh::FEMContext::FENeeded::value_getter
void(FEMContext::* value_getter)(unsigned int, value_base *&, unsigned short) const
Definition: fem_context.h:1087
libMesh::FEMContext::_custom_solution
const NumericVector< Number > * _custom_solution
Data with which to do algebra reinitialization.
Definition: fem_context.h:1045
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::FEMContext::_mesh_x_var
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: fem_context.h:1005
libMesh::FEMContext::attach_quadrature_rules
void attach_quadrature_rules()
Helper function for attaching quadrature rules.
Definition: fem_context.C:112
libMesh::DenseSubVector< Number >
libMesh::FEMContext::use_default_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.
Definition: fem_context.C:134
libMesh::FEMContext::_update_time_from_system
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...
Definition: fem_context.C:1864
libMesh::FEMContext::set_mesh_system
virtual void set_mesh_system(System *sys)
Tells the FEMContext that system sys contains the isoparametric Lagrangian variables which correspond...
Definition: fem_context.h:830
libMesh::FEMContext::diff_subsolution_getter
const typedef DenseSubVector< Number > &(DiffContext::* diff_subsolution_getter)(unsigned int) const
Helper typedef to simplify refactoring.
Definition: fem_context.h:1071
libMesh::FEMContext::FENeeded::value_base
FEGenericBase< value_shape > value_base
Definition: fem_context.h:1086
libMesh::FEMContext::NONE
Definition: fem_context.h:958
libMesh::FEMContext::point_hessian
Tensor point_hessian(unsigned int var, const Point &p) const
Definition: fem_context.C:931
libMesh::FEMContext::elem_fe_reinit
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
Reinitializes interior FE objects on the current geometric element.
Definition: fem_context.C:1436
libMesh::FEMContext::point_accel
void point_accel(unsigned int var, const Point &p, OutputType &u) const
libMesh::FEMContext::edge_fe_reinit
virtual void edge_fe_reinit()
Reinitializes edge FE objects on the current geometric element.
Definition: fem_context.C:1473
libMesh::FEMContext::get_element_fe
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.
Definition: fem_context.h:275
libMesh::FEMContext::_side_qrule
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 ...
Definition: fem_context.h:1192
libMesh::FEMContext::interior_gradient
Gradient interior_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:426
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::FEMContext::set_mesh_y_var
void set_mesh_y_var(unsigned int y_var)
Accessor for y-variable of moving mesh System.
Definition: fem_context.h:870
libMesh::FEMContext::DOFS_ONLY
Definition: fem_context.h:959
libMesh::FEMContext::set_algebraic_type
void set_algebraic_type(const AlgebraicType atype)
Setting which determines whether to initialize algebraic structures (elem_*) on each element and set ...
Definition: fem_context.h:973
libMesh::FEMContext::FENeeded::value_shape
TensorTools::MakeReal< OutputType >::type value_shape
Definition: fem_context.h:1085
libMesh::FEMContext::FENeeded::hess_base
FEGenericBase< hess_shape > hess_base
Definition: fem_context.h:1096
libMesh::FEMContext::has_side_boundary_id
bool has_side_boundary_id(boundary_id_type id) const
Reports if the boundary id is found on the current side.
Definition: fem_context.C:254
libMesh::FEMContext::side_value
Number side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:620
libMesh::FEMContext::algebraic_type
AlgebraicType algebraic_type() const
Definition: fem_context.h:979
libMesh::FEMContext::get_edge_qrule
const QBase & get_edge_qrule() const
Accessor for element edge quadrature rule.
Definition: fem_context.h:819
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::FEMContext::get_side_qrule
const QBase & get_side_qrule() const
Accessor for element side quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:797
libMesh::FEMContext::get_dim
unsigned char get_dim() const
Accessor for cached mesh dimension.
Definition: fem_context.h:924
libMesh::FEMContext::_mesh_z_var
unsigned int _mesh_z_var
Definition: fem_context.h:1005
libMesh::FEMContext::AlgebraicType
AlgebraicType
Enum describing what data to use when initializing algebraic structures on each element.
Definition: fem_context.h:958
libMesh::FEMContext::fixed_interior_value
Number fixed_interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1018
libMesh::FEMContext::find_hardest_fe_type
FEType find_hardest_fe_type()
Helper function for creating quadrature rules.
Definition: fem_context.C:80
libMesh::FEMContext::get_side_fe
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...
Definition: fem_context.h:312
libMesh::FEMContext::some_value
void some_value(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_value methods.
Definition: fem_context.C:279
libMesh::FEMContext::elem_dimensions
const std::set< unsigned char > & elem_dimensions() const
Definition: fem_context.h:938
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEMContext::OLD
Definition: fem_context.h:962
libMesh::FEMContext::get_mesh_z_var
unsigned int get_mesh_z_var() const
Accessor for z-variable of moving mesh System.
Definition: fem_context.h:876
libMesh::FEMContext::_real_grad_fe_is_inf
bool _real_grad_fe_is_inf
Definition: fem_context.h:1052
libMesh::FEMContext::elem_edge_reinit
virtual void elem_edge_reinit(Real theta) override
Resets the current time in the context.
Definition: fem_context.C:1411
libMesh::FEMContext::get_edge
unsigned char get_edge() const
Accessor for current edge of Elem object.
Definition: fem_context.h:916
libMesh::FEMContext::interior_div
void interior_div(unsigned int var, unsigned int qp, OutputType &div_u) const
Definition: fem_context.C:588
libMesh::FEMContext::fixed_interior_hessian
Tensor fixed_interior_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1068
libMesh::FEMContext::_element_fe_var
std::vector< std::vector< FEAbstract * > > _element_fe_var
Pointers to the same finite element objects, but indexed by variable number.
Definition: fem_context.h:1147
libMesh::FEMContext::CURRENT
Definition: fem_context.h:961
libMesh::FEMContext::interior_value
Number interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:371
libMesh::FEMContext::set_mesh_z_var
void set_mesh_z_var(unsigned int z_var)
Accessor for z-variable of moving mesh System.
Definition: fem_context.h:884
libMesh::FEMContext::fixed_point_hessian
Tensor fixed_point_hessian(unsigned int var, const Point &p) const
Definition: fem_context.C:1265
int
void ErrorVector unsigned int
Definition: adjoints_ex3.C:360
libMesh::FEMContext::_dim
unsigned char _dim
Cached dimension of largest dimension element in this mesh.
Definition: fem_context.h:1165
libMesh::FEMContext::set_mesh_x_var
void set_mesh_x_var(unsigned int x_var)
Accessor for x-variable of moving mesh System.
Definition: fem_context.h:856
libMesh::FEMContext::_edge_fe_var
std::vector< FEAbstract * > _edge_fe_var
Definition: fem_context.h:1149
libMesh::FEMContext::elem_reinit
virtual void elem_reinit(Real theta) override
Resets the current time in the context.
Definition: fem_context.C:1372
libMesh::FEMContext::side_hessians
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...
Definition: fem_context.C:790
libMesh::FEMContext::fixed_interior_gradient
Gradient fixed_interior_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1041
libMesh::FEMContext::get_elem
Elem & get_elem()
Accessor for current Elem object.
Definition: fem_context.h:903
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62