fem_context.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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 
37 namespace libMesh
38 {
39 
40 // Forward Declarations
41 class BoundaryInfo;
42 class Elem;
43 template <typename T> class FEGenericBase;
44 typedef FEGenericBase<Real> FEBase;
45 class QBase;
46 class Point;
47 template <typename T> class NumericVector;
48 
61 class FEMContext : public DiffContext
62 {
63 public:
64 
68  explicit
69  FEMContext (const System & sys);
70 
75  explicit
76  FEMContext (const System & sys, int extra_quadrature_order);
77 
81  virtual ~FEMContext ();
82 
87 
94 #ifdef LIBMESH_ENABLE_DEPRECATED
95  std::vector<boundary_id_type> side_boundary_ids() const;
96 #endif
97 
101  void side_boundary_ids(std::vector<boundary_id_type> & vec_to_fill) const;
102 
109  Number interior_value(unsigned int var, unsigned int qp) const;
110 
117  Number side_value(unsigned int var, unsigned int qp) const;
118 
125  Number point_value(unsigned int var, const Point & p) const;
126 
133  Gradient interior_gradient(unsigned int var, unsigned int qp) const;
134 
141  Gradient side_gradient(unsigned int var, unsigned int qp) const;
142 
149  Gradient point_gradient(unsigned int var, const Point & p) const;
150 
151 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
152 
158  Tensor interior_hessian(unsigned int var, unsigned int qp) const;
159 
166  Tensor side_hessian(unsigned int var, unsigned int qp) const;
167 
174  Tensor point_hessian(unsigned int var, const Point & p) const;
175 
176 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
177 
184  Number fixed_interior_value(unsigned int var, unsigned int qp) const;
185 
192  Number fixed_side_value(unsigned int var, unsigned int qp) const;
193 
200  Number fixed_point_value(unsigned int var, const Point & p) const;
201 
208  Gradient fixed_interior_gradient(unsigned int var, unsigned int qp) const;
209 
216  Gradient fixed_side_gradient(unsigned int var, unsigned int qp) const;
217 
224  Gradient fixed_point_gradient(unsigned int var, const Point & p) const;
225 
226 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
227 
233  Tensor fixed_interior_hessian(unsigned int var, unsigned int qp) const;
234 
241  Tensor fixed_side_hessian(unsigned int var, unsigned int qp) const;
242 
249  Tensor fixed_point_hessian (unsigned int var, const Point & p) const;
250 
251 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
252 
261  template<typename OutputShape>
262  void get_element_fe( unsigned int var, FEGenericBase<OutputShape> *& fe ) const
263  { this->get_element_fe<OutputShape>(var,fe,this->get_dim()); }
264 
273  FEBase * get_element_fe( unsigned int var ) const
274  { return this->get_element_fe(var,this->get_dim()); }
275 
280  template<typename OutputShape>
281  void get_element_fe( unsigned int var, FEGenericBase<OutputShape> *& fe,
282  unsigned short dim ) const;
283 
288  FEBase * get_element_fe( unsigned int var, unsigned short dim ) const;
289 
298  template<typename OutputShape>
299  void get_side_fe( unsigned int var, FEGenericBase<OutputShape> *& fe ) const
300  { this->get_side_fe<OutputShape>(var,fe,this->get_dim()); }
301 
310  FEBase * get_side_fe( unsigned int var ) const
311  { return this->get_side_fe(var,this->get_dim()); }
312 
317  template<typename OutputShape>
318  void get_side_fe( unsigned int var, FEGenericBase<OutputShape> *& fe,
319  unsigned short dim ) const;
320 
325  FEBase * get_side_fe( unsigned int var, unsigned short dim ) const;
326 
330  template<typename OutputShape>
331  void get_edge_fe( unsigned int var, FEGenericBase<OutputShape> *& fe ) const;
332 
336  FEBase * get_edge_fe( unsigned int var ) const;
337 
344  template<typename OutputType>
345  void interior_value(unsigned int var,
346  unsigned int qp,
347  OutputType & u) const;
348 
353  template<typename OutputType>
354  void interior_values(unsigned int var,
355  const NumericVector<Number> & _system_vector,
356  std::vector<OutputType> & interior_values_vector) const;
357 
364  template<typename OutputType>
365  void side_value(unsigned int var,
366  unsigned int qp,
367  OutputType & u) const;
368 
373  template<typename OutputType>
374  void side_values(unsigned int var,
375  const NumericVector<Number> & _system_vector,
376  std::vector<OutputType> & side_values_vector) const;
377 
387  template<typename OutputType>
388  void point_value(unsigned int var,
389  const Point & p,
390  OutputType & u,
391  const Real tolerance = TOLERANCE) const;
392 
399  template<typename OutputType>
400  void interior_gradient(unsigned int var, unsigned int qp,
401  OutputType & du) const;
402 
409  template<typename OutputType>
410  void interior_gradients(unsigned int var,
411  const NumericVector<Number> & _system_vector,
412  std::vector<OutputType> & interior_gradients_vector) const;
413 
420  template<typename OutputType>
421  void side_gradient(unsigned int var,
422  unsigned int qp,
423  OutputType & du) const;
424 
431  template<typename OutputType>
432  void side_gradients(unsigned int var,
433  const NumericVector<Number> & _system_vector,
434  std::vector<OutputType> & side_gradients_vector) const;
435 
445  template<typename OutputType>
446  void point_gradient(unsigned int var,
447  const Point & p,
448  OutputType & grad_u,
449  const Real tolerance = TOLERANCE) const;
450 
451 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
452 
458  template<typename OutputType>
459  void interior_hessian(unsigned int var,
460  unsigned int qp,
461  OutputType & d2u) const;
462 
468  template<typename OutputType>
469  void interior_hessians(unsigned int var,
470  const NumericVector<Number> & _system_vector,
471  std::vector<OutputType> & d2u_vals) const;
472 
479  template<typename OutputType>
480  void side_hessian(unsigned int var,
481  unsigned int qp,
482  OutputType & d2u) const;
483 
489  template<typename OutputType>
490  void side_hessians(unsigned int var,
491  const NumericVector<Number> & _system_vector,
492  std::vector<OutputType> & d2u_vals) const;
493 
503  template<typename OutputType>
504  void point_hessian(unsigned int var,
505  const Point & p,
506  OutputType & hess_u,
507  const Real tolerance = TOLERANCE) const;
508 
509 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
510 
516  template<typename OutputType>
517  void interior_rate(unsigned int var,
518  unsigned int qp,
519  OutputType & u) const;
520 
525  template<typename OutputType>
526  void side_rate(unsigned int var,
527  unsigned int qp,
528  OutputType & u) const;
529 
534  template<typename OutputType>
535  void point_rate(unsigned int var,
536  const Point & p,
537  OutputType & u) const;
538 
544  template<typename OutputType>
545  void interior_accel(unsigned int var,
546  unsigned int qp,
547  OutputType & u) const;
548 
553  template<typename OutputType>
554  void side_accel(unsigned int var,
555  unsigned int qp,
556  OutputType & u) const;
557 
562  template<typename OutputType>
563  void point_accel(unsigned int var,
564  const Point & p,
565  OutputType & u) const;
566 
573  template<typename OutputType>
574  void fixed_interior_value(unsigned int var,
575  unsigned int qp,
576  OutputType & u) const;
577 
584  template<typename OutputType>
585  void fixed_side_value(unsigned int var,
586  unsigned int qp,
587  OutputType & u) const;
588 
598  template<typename OutputType>
599  void fixed_point_value(unsigned int var,
600  const Point & p,
601  OutputType & u,
602  const Real tolerance = TOLERANCE) const;
603 
610  template<typename OutputType>
611  void fixed_interior_gradient(unsigned int var,
612  unsigned int qp,
613  OutputType & grad_u) const;
614 
621  template<typename OutputType>
622  void fixed_side_gradient(unsigned int var,
623  unsigned int qp,
624  OutputType & grad_u) const;
625 
635  template<typename OutputType>
636  void fixed_point_gradient(unsigned int var,
637  const Point & p,
638  OutputType & grad_u,
639  const Real tolerance = TOLERANCE) const;
640 
641 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
642 
648  template<typename OutputType>
649  void fixed_interior_hessian(unsigned int var,
650  unsigned int qp,
651  OutputType & hess_u) const;
652 
659  template<typename OutputType>
660  void fixed_side_hessian(unsigned int var,
661  unsigned int qp,
662  OutputType & hess_u) const;
663 
673  template<typename OutputType>
674  void fixed_point_hessian(unsigned int var,
675  const Point & p,
676  OutputType & hess_u,
677  const Real tolerance = TOLERANCE) const;
678 
679 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
680 
685  template<typename OutputType>
686  void interior_curl(unsigned int var,
687  unsigned int qp,
688  OutputType & curl_u) const;
689 
697  template<typename OutputType>
698  void point_curl(unsigned int var,
699  const Point & p,
700  OutputType & curl_u,
701  const Real tolerance = TOLERANCE) const;
702 
707  template<typename OutputType>
708  void interior_div(unsigned int var,
709  unsigned int qp,
710  OutputType & div_u) const;
711 
712  // should be protected:
718  virtual void elem_reinit(Real theta) override;
719 
725  virtual void elem_side_reinit(Real theta) override;
726 
732  virtual void elem_edge_reinit(Real theta) override;
733 
738  virtual void nonlocal_reinit(Real theta) override;
739 
743  virtual void pre_fe_reinit(const System &,
744  const Elem * e);
745 
749  virtual void elem_fe_reinit(const std::vector<Point> * const pts = nullptr);
750 
754  virtual void side_fe_reinit();
755 
759  virtual void edge_fe_reinit();
760 
765  const QBase & get_element_qrule() const
766  { return this->get_element_qrule(this->get_elem_dim()); }
767 
772  const QBase & get_side_qrule() const
773  { return this->get_side_qrule(this->get_elem_dim()); }
774 
778  const QBase & get_element_qrule( unsigned short dim ) const
779  { libmesh_assert(_element_qrule[dim]);
780  return *(this->_element_qrule[dim]); }
781 
785  const QBase & get_side_qrule( unsigned short dim ) const
786  {
787  libmesh_assert(_side_qrule[dim]);
788  return *(this->_side_qrule[dim]);
789  }
790 
794  const QBase & get_edge_qrule() const
795  { return *(this->_edge_qrule); }
796 
805  virtual void set_mesh_system(System * sys)
806  { this->_mesh_sys = sys; }
807 
811  const System * get_mesh_system() const
812  { return this->_mesh_sys; }
813 
818  { return this->_mesh_sys; }
819 
823  unsigned int get_mesh_x_var() const
824  { return _mesh_x_var; }
825 
831  void set_mesh_x_var(unsigned int x_var)
832  { _mesh_x_var = x_var; }
833 
837  unsigned int get_mesh_y_var() const
838  { return _mesh_y_var; }
839 
845  void set_mesh_y_var(unsigned int y_var)
846  { _mesh_y_var = y_var; }
847 
851  unsigned int get_mesh_z_var() const
852  { return _mesh_z_var; }
853 
859  void set_mesh_z_var(unsigned int z_var)
860  { _mesh_z_var = z_var; }
861 
865  bool has_elem() const
866  { return (this->_elem != nullptr); }
867 
871  const Elem & get_elem() const
872  { libmesh_assert(this->_elem);
873  return *(this->_elem); }
874 
879  { libmesh_assert(this->_elem);
880  return *(const_cast<Elem *>(this->_elem)); }
881 
885  unsigned char get_side() const
886  { return side; }
887 
891  unsigned char get_edge() const
892  { return edge; }
893 
899  unsigned char get_dim() const
900  { return this->_dim; }
901 
906  unsigned char get_elem_dim() const
907  { return _elem_dim; }
908 
913  const std::set<unsigned char> & elem_dimensions() const
914  { return _elem_dims; }
915 
921  void elem_position_set(Real theta);
922 
927  void elem_position_get();
928 
933  enum AlgebraicType { NONE = 0, // Do not reinitialize dof_indices
934  DOFS_ONLY, // Reinitialize dof_indices, not
935  // algebraic structures
936  CURRENT, // Use dof_indices, current solution
937  OLD, // Use old_dof_indices, custom solution
938  OLD_DOFS_ONLY}; // Reinitialize old_dof_indices, not
939  // algebraic structures
940 
949  { _atype = atype; }
950 
951  /*
952  * Get the current AlgebraicType setting
953  */
954  AlgebraicType algebraic_type() const { return _atype; }
955 
964  { _custom_solution = custom_sol; }
965 
970 
975 
979  unsigned char side;
980 
984  unsigned char edge;
985 
989  template<typename OutputShape>
991  const Point & p,
992  const Real tolerance = TOLERANCE) const;
993 
994 protected:
995 
1000 
1005 
1006  mutable std::unique_ptr<FEGenericBase<Real>> _real_fe;
1007  mutable std::unique_ptr<FEGenericBase<RealGradient>> _real_grad_fe;
1008 
1009 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1010  mutable bool _real_fe_is_inf;
1011  mutable bool _real_grad_fe_is_inf;
1012 #endif
1013 
1014  template<typename OutputShape>
1015  FEGenericBase<OutputShape> * cached_fe( const unsigned int elem_dim,
1016  const FEType fe_type ) const;
1017 
1021  void set_elem( const Elem * e );
1022 
1023  // gcc-3.4, oracle 12.3 require this typedef to be public
1024  // in order to use it in a return type
1025 public:
1026 
1030  typedef const DenseSubVector<Number> & (DiffContext::*diff_subsolution_getter)(unsigned int) const;
1031 
1032 protected:
1036  template <typename OutputType>
1037  struct FENeeded
1038  {
1039  // Rank decrementer helper types
1042 
1043  // Typedefs for "Value getter" function pointer
1046  typedef void (FEMContext::*value_getter) (unsigned int, value_base *&, unsigned short) const;
1047 
1048  // Typedefs for "Grad getter" function pointer
1051  typedef void (FEMContext::*grad_getter) (unsigned int, grad_base *&, unsigned short) const;
1052 
1053  // Typedefs for "Hessian getter" function pointer
1056  typedef void (FEMContext::*hess_getter) (unsigned int, hess_base *&, unsigned short) const;
1057  };
1058 
1059 
1060 
1065  template<typename OutputType,
1066  typename FENeeded<OutputType>::value_getter fe_getter,
1067  diff_subsolution_getter subsolution_getter>
1068  void some_value(unsigned int var, unsigned int qp, OutputType & u) const;
1069 
1074  template<typename OutputType,
1075  typename FENeeded<OutputType>::grad_getter fe_getter,
1076  diff_subsolution_getter subsolution_getter>
1077  void some_gradient(unsigned int var, unsigned int qp, OutputType & u) const;
1078 
1079 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1080 
1084  template<typename OutputType,
1085  typename FENeeded<OutputType>::hess_getter fe_getter,
1086  diff_subsolution_getter subsolution_getter>
1087  void some_hessian(unsigned int var, unsigned int qp, OutputType & u) const;
1088 #endif
1089 
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;
1098 
1099 
1106  std::vector<std::vector<FEAbstract *>> _element_fe_var;
1107  std::vector<std::vector<FEAbstract *>> _side_fe_var;
1108  std::vector<FEAbstract *> _edge_fe_var;
1109 
1115 
1119  const Elem * _elem;
1120 
1124  unsigned char _dim;
1125 
1129  unsigned char _elem_dim;
1130 
1135  std::set<unsigned char> _elem_dims;
1136 
1143  std::vector<std::unique_ptr<QBase>> _element_qrule;
1144 
1151  std::vector<std::unique_ptr<QBase>> _side_qrule;
1152 
1160  std::unique_ptr<QBase> _edge_qrule;
1161 
1166 
1167 private:
1171  void init_internal_data(const System & sys);
1172 
1181  void _do_elem_position_set(Real theta);
1182 
1188  void _update_time_from_system(Real theta);
1189 };
1190 
1191 
1192 
1193 // ------------------------------------------------------------
1194 // FEMContext inline methods
1195 
1196 inline
1198 {
1199  if (_mesh_sys)
1200  this->_do_elem_position_set(theta);
1201 }
1202 
1203 template<typename OutputShape>
1204 inline
1206  unsigned short dim ) const
1207 {
1208  libmesh_assert( !_element_fe_var[dim].empty() );
1209  libmesh_assert_less ( var, (_element_fe_var[dim].size() ) );
1210  fe = cast_ptr<FEGenericBase<OutputShape> *>( (_element_fe_var[dim][var] ) );
1211 }
1212 
1213 inline
1214 FEBase * FEMContext::get_element_fe( unsigned int var, unsigned short dim ) const
1215 {
1216  libmesh_assert( !_element_fe_var[dim].empty() );
1217  libmesh_assert_less ( var, (_element_fe_var[dim].size() ) );
1218  return cast_ptr<FEBase *>( (_element_fe_var[dim][var] ) );
1219 }
1220 
1221 template<typename OutputShape>
1222 inline
1224  unsigned short dim ) const
1225 {
1226  libmesh_assert( !_side_fe_var[dim].empty() );
1227  libmesh_assert_less ( var, (_side_fe_var[dim].size() ) );
1228  fe = cast_ptr<FEGenericBase<OutputShape> *>( (_side_fe_var[dim][var] ) );
1229 }
1230 
1231 inline
1232 FEBase * FEMContext::get_side_fe( unsigned int var, unsigned short dim ) const
1233 {
1234  libmesh_assert( !_side_fe_var[dim].empty() );
1235  libmesh_assert_less ( var, (_side_fe_var[dim].size() ) );
1236  return cast_ptr<FEBase *>( (_side_fe_var[dim][var] ) );
1237 }
1238 
1239 template<typename OutputShape>
1240 inline
1241 void FEMContext::get_edge_fe( unsigned int var, FEGenericBase<OutputShape> *& fe ) const
1242 {
1243  libmesh_assert_less ( var, _edge_fe_var.size() );
1244  fe = cast_ptr<FEGenericBase<OutputShape> *>( _edge_fe_var[var] );
1245 }
1246 
1247 inline
1248 FEBase * FEMContext::get_edge_fe( unsigned int var ) const
1249 {
1250  libmesh_assert_less ( var, _edge_fe_var.size() );
1251  return cast_ptr<FEBase *>( _edge_fe_var[var] );
1252 }
1253 
1254 
1255 } // namespace libMesh
1256 
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.
Definition: fe_type.h:179
unsigned char get_edge() const
Definition: fem_context.h:891
virtual void nonlocal_reinit(Real theta) override
Definition: fem_context.C:1352
void(FEMContext::* value_getter)(unsigned int, value_base *&, unsigned short) const
Definition: fem_context.h:1046
Number fixed_side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1029
Number side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:555
Tensor fixed_point_hessian(unsigned int var, const Point &p) const
Definition: fem_context.C:1200
void set_mesh_z_var(unsigned int z_var)
Definition: fem_context.h:859
void point_rate(unsigned int var, const Point &p, OutputType &u) const
Number interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:306
void interior_div(unsigned int var, unsigned int qp, OutputType &div_u) const
Definition: fem_context.C:523
void elem_position_set(Real theta)
Definition: fem_context.h:1197
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Definition: fem_context.h:299
virtual void pre_fe_reinit(const System &, const Elem *e)
Definition: fem_context.C:1552
unsigned int get_mesh_x_var() const
Definition: fem_context.h:823
void set_mesh_x_var(unsigned int x_var)
Definition: fem_context.h:831
const NumericVector< Number > * _custom_solution
Definition: fem_context.h:1004
void set_elem(const Elem *e)
Definition: fem_context.C:1765
void get_edge_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Definition: fem_context.h:1241
unsigned int _mesh_z_var
Definition: fem_context.h:974
std::unique_ptr< FEGenericBase< RealGradient > > _real_grad_fe
Definition: fem_context.h:1007
unsigned char get_elem_dim() const
Definition: fem_context.h:906
const Elem & get_elem() const
Definition: fem_context.h:871
TensorTools::MakeReal< OutputType >::type value_shape
Definition: fem_context.h:1044
const QBase & get_element_qrule(unsigned short dim) const
Definition: fem_context.h:778
const DenseSubVector< Number > &(DiffContext::* diff_subsolution_getter)(unsigned int) const
Definition: fem_context.h:1030
FEGenericBase< value_shape > value_base
Definition: fem_context.h:1045
virtual void elem_edge_reinit(Real theta) override
Definition: fem_context.C:1337
The base class for all geometric element types.
Definition: elem.h:100
void set_algebraic_type(const AlgebraicType atype)
Definition: fem_context.h:948
static const Real TOLERANCE
FEGenericBase< hess_shape > hess_base
Definition: fem_context.h:1055
void interior_gradients(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_gradients_vector) const
Definition: fem_context.C:388
bool has_side_boundary_id(boundary_id_type id) const
Definition: fem_context.C:189
virtual ~FEMContext()
Definition: fem_context.C:183
unsigned char get_side() const
Definition: fem_context.h:885
Number point_value(unsigned int var, const Point &p) const
Definition: fem_context.C:768
DiffContext(const System &)
Definition: diff_context.C:30
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _element_fe
Definition: fem_context.h:1095
unsigned char get_dim() const
Definition: fem_context.h:899
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
Definition: fem_context.C:1362
std::unique_ptr< FEGenericBase< Real > > _real_fe
Definition: fem_context.h:1006
Number fixed_interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:953
virtual void set_mesh_system(System *sys)
Definition: fem_context.h:805
Tensor fixed_side_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1078
unsigned char edge
Definition: fem_context.h:984
std::set< unsigned char > _elem_dims
Definition: fem_context.h:1135
Gradient interior_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:361
void point_curl(unsigned int var, const Point &p, OutputType &curl_u, const Real tolerance=TOLERANCE) const
Definition: fem_context.C:917
int8_t boundary_id_type
Definition: id_types.h:51
void side_values(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_values_vector) const
Definition: fem_context.C:578
unsigned char _elem_dim
Definition: fem_context.h:1129
Tensor side_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:696
virtual void side_fe_reinit()
Definition: fem_context.C:1383
void interior_hessians(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
Definition: fem_context.C:450
std::map< FEType, std::unique_ptr< FEAbstract > > _edge_fe
Definition: fem_context.h:1097
void _update_time_from_system(Real theta)
Definition: fem_context.C:1774
void init_internal_data(const System &sys)
Definition: fem_context.C:79
Gradient fixed_point_gradient(unsigned int var, const Point &p) const
Definition: fem_context.C:1149
const QBase & get_side_qrule(unsigned short dim) const
Definition: fem_context.h:785
void interior_values(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_values_vector) const
Definition: fem_context.C:326
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
void some_value(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:214
void some_hessian(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:276
Tensor interior_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:426
Tensor point_hessian(unsigned int var, const Point &p) const
Definition: fem_context.C:866
Gradient fixed_interior_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:976
unsigned int _mesh_x_var
Definition: fem_context.h:974
Used by the Mesh to keep track of boundary nodes and elements.
Definition: boundary_info.h:57
TensorTools::MakeReal< Rank1Decrement >::type grad_shape
Definition: fem_context.h:1049
TensorTools::MakeReal< Rank2Decrement >::type hess_shape
Definition: fem_context.h:1054
AlgebraicType _atype
Definition: fem_context.h:999
const BoundaryInfo & _boundary_info
Definition: fem_context.h:1114
std::vector< std::vector< FEAbstract * > > _element_fe_var
Definition: fem_context.h:1106
FEGenericBase< grad_shape > grad_base
Definition: fem_context.h:1050
void(FEMContext::* hess_getter)(unsigned int, hess_base *&, unsigned short) const
Definition: fem_context.h:1056
FEBase * get_side_fe(unsigned int var) const
Definition: fem_context.h:310
std::vector< boundary_id_type > side_boundary_ids() const
Definition: fem_context.C:196
FEGenericBase< Real > FEBase
TensorTools::DecrementRank< Rank1Decrement >::type Rank2Decrement
Definition: fem_context.h:1041
virtual void elem_side_reinit(Real theta) override
Definition: fem_context.C:1322
Number fixed_point_value(unsigned int var, const Point &p) const
Definition: fem_context.C:1103
unsigned char _dim
Definition: fem_context.h:1124
virtual void edge_fe_reinit()
Definition: fem_context.C:1399
std::vector< FEAbstract * > _edge_fe_var
Definition: fem_context.h:1108
void interior_rate(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1253
std::vector< std::unique_ptr< QBase > > _side_qrule
Definition: fem_context.h:1151
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _side_fe
Definition: fem_context.h:1096
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< std::unique_ptr< QBase > > _element_qrule
Definition: fem_context.h:1143
void side_accel(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1287
unsigned char side
Definition: fem_context.h:979
void interior_curl(unsigned int var, unsigned int qp, OutputType &curl_u) const
Definition: fem_context.C:492
TensorTools::DecrementRank< OutputType >::type Rank1Decrement
Definition: fem_context.h:1040
void set_mesh_y_var(unsigned int y_var)
Definition: fem_context.h:845
std::unique_ptr< QBase > _edge_qrule
Definition: fem_context.h:1160
FEMContext(const System &sys)
Definition: fem_context.C:37
const System * get_mesh_system() const
Definition: fem_context.h:811
void side_gradients(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_gradients_vector) const
Definition: fem_context.C:658
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Definition: fem_context.h:262
unsigned int _mesh_y_var
Definition: fem_context.h:974
Gradient fixed_side_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1053
unsigned int get_mesh_y_var() const
Definition: fem_context.h:837
void(FEMContext::* grad_getter)(unsigned int, grad_base *&, unsigned short) const
Definition: fem_context.h:1051
AlgebraicType algebraic_type() const
Definition: fem_context.h:954
const QBase & get_edge_qrule() const
Definition: fem_context.h:794
void some_gradient(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:243
void set_custom_solution(const NumericVector< Number > *custom_sol)
Definition: fem_context.h:963
std::vector< std::vector< FEAbstract * > > _side_fe_var
Definition: fem_context.h:1107
FEGenericBase< OutputShape > * build_new_fe(const FEGenericBase< OutputShape > *fe, const Point &p, const Real tolerance=TOLERANCE) const
Definition: fem_context.C:1862
Gradient point_gradient(unsigned int var, const Point &p) const
Definition: fem_context.C:814
const Elem * _elem
Definition: fem_context.h:1119
FEGenericBase< OutputShape > * cached_fe(const unsigned int elem_dim, const FEType fe_type) const
void _do_elem_position_set(Real theta)
Definition: fem_context.C:1476
bool has_elem() const
Definition: fem_context.h:865
Gradient side_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:613
FEBase * get_element_fe(unsigned int var) const
Definition: fem_context.h:273
void side_hessians(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
Definition: fem_context.C:725
void side_rate(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1265
A geometric point in (x,y,z) space.
Definition: point.h:38
System * get_mesh_system()
Definition: fem_context.h:817
Base class for all quadrature families and orders.
Definition: quadrature.h:62
const QBase & get_element_qrule() const
Definition: fem_context.h:765
virtual void elem_reinit(Real theta) override
Definition: fem_context.C:1298
unsigned int get_mesh_z_var() const
Definition: fem_context.h:851
const std::set< unsigned char > & elem_dimensions() const
Definition: fem_context.h:913
void interior_accel(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1275
const QBase & get_side_qrule() const
Definition: fem_context.h:772
Tensor fixed_interior_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1003