libMesh
fe_abstract.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_FE_ABSTRACT_H
21 #define LIBMESH_FE_ABSTRACT_H
22 
23 // Local includes
24 #include "libmesh/reference_counted_object.h"
25 #include "libmesh/point.h"
26 #include "libmesh/vector_value.h"
27 #include "libmesh/fe_type.h"
28 #include "libmesh/fe_map.h"
29 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
30 #include "libmesh/tensor_value.h"
31 #endif
32 
33 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS
34 namespace libMesh
35 {
36 enum ElemType : int;
37 }
38 #else
39 #include "libmesh/enum_elem_type.h"
40 #endif
41 
42 // C++ includes
43 #include <cstddef>
44 #include <vector>
45 #include <memory>
46 
47 namespace libMesh
48 {
49 
50 
51 // forward declarations
52 template <typename T> class DenseMatrix;
53 template <typename T> class DenseVector;
54 class BoundaryInfo;
55 class DofConstraints;
56 class DofMap;
57 class Elem;
58 class MeshBase;
59 template <typename T> class NumericVector;
60 class QBase;
61 
62 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
63 class NodeConstraints;
64 #endif
65 
66 #ifdef LIBMESH_ENABLE_PERIODIC
67 class PeriodicBoundaries;
68 class PointLocatorBase;
69 #endif
70 
71 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
72 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
73 class InfFE;
74 #endif
75 
76 
77 
100 class FEAbstract : public ReferenceCountedObject<FEAbstract>
101 {
102 protected:
103 
109  FEAbstract (const unsigned int dim,
110  const FEType & fet);
111 
112 public:
113 
117  virtual ~FEAbstract();
118 
125  static std::unique_ptr<FEAbstract> build (const unsigned int dim,
126  const FEType & type);
127 
142  virtual void reinit (const Elem * elem,
143  const std::vector<Point> * const pts = nullptr,
144  const std::vector<Real> * const weights = nullptr) = 0;
145 
154  virtual void reinit (const Elem * elem,
155  const unsigned int side,
156  const Real tolerance = TOLERANCE,
157  const std::vector<Point> * const pts = nullptr,
158  const std::vector<Real> * const weights = nullptr) = 0;
159 
168  virtual void edge_reinit (const Elem * elem,
169  const unsigned int edge,
170  const Real tolerance = TOLERANCE,
171  const std::vector<Point> * pts = nullptr,
172  const std::vector<Real> * weights = nullptr) = 0;
173 
178  virtual void side_map (const Elem * elem,
179  const Elem * side,
180  const unsigned int s,
181  const std::vector<Point> & reference_side_points,
182  std::vector<Point> & reference_points) = 0;
183 
191  static bool on_reference_element(const Point & p,
192  const ElemType t,
193  const Real eps = TOLERANCE);
198  static void get_refspace_nodes(const ElemType t,
199  std::vector<Point> & nodes);
200 
201 
202 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
203 
207  static void compute_node_constraints (NodeConstraints & constraints,
208  const Elem * elem);
209 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
210 
211 #ifdef LIBMESH_ENABLE_PERIODIC
212 
213 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
214 
218  static void compute_periodic_node_constraints (NodeConstraints & constraints,
219  const PeriodicBoundaries & boundaries,
220  const MeshBase & mesh,
221  const PointLocatorBase * point_locator,
222  const Elem * elem);
223 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
224 
225 #endif // LIBMESH_ENABLE_PERIODIC
226 
230  unsigned int get_dim() const
231  { return dim; }
232 
237  const std::vector<Point> & get_xyz() const
238  { calculate_map = true; return this->_fe_map->get_xyz(); }
239 
244  const std::vector<Real> & get_JxW() const
245  { calculate_map = true; return this->_fe_map->get_JxW(); }
246 
251  const std::vector<RealGradient> & get_dxyzdxi() const
252  { calculate_map = true; return this->_fe_map->get_dxyzdxi(); }
253 
258  const std::vector<RealGradient> & get_dxyzdeta() const
259  { calculate_map = true; return this->_fe_map->get_dxyzdeta(); }
260 
265  const std::vector<RealGradient> & get_dxyzdzeta() const
266  { return _fe_map->get_dxyzdzeta(); }
267 
268 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
269 
273  const std::vector<RealGradient> & get_d2xyzdxi2() const
274  { calculate_map = true; return this->_fe_map->get_d2xyzdxi2(); }
275 
279  const std::vector<RealGradient> & get_d2xyzdeta2() const
280  { calculate_map = true; return this->_fe_map->get_d2xyzdeta2(); }
281 
285  const std::vector<RealGradient> & get_d2xyzdzeta2() const
286  { calculate_map = true; return this->_fe_map->get_d2xyzdzeta2(); }
287 
291  const std::vector<RealGradient> & get_d2xyzdxideta() const
292  { calculate_map = true; return this->_fe_map->get_d2xyzdxideta(); }
293 
297  const std::vector<RealGradient> & get_d2xyzdxidzeta() const
298  { calculate_map = true; return this->_fe_map->get_d2xyzdxidzeta(); }
299 
303  const std::vector<RealGradient> & get_d2xyzdetadzeta() const
304  { calculate_map = true; return this->_fe_map->get_d2xyzdetadzeta(); }
305 
306 #endif
307 
312  const std::vector<Real> & get_dxidx() const
313  { calculate_map = true; return this->_fe_map->get_dxidx(); }
314 
319  const std::vector<Real> & get_dxidy() const
320  { calculate_map = true; return this->_fe_map->get_dxidy(); }
321 
326  const std::vector<Real> & get_dxidz() const
327  { calculate_map = true; return this->_fe_map->get_dxidz(); }
328 
333  const std::vector<Real> & get_detadx() const
334  { calculate_map = true; return this->_fe_map->get_detadx(); }
335 
340  const std::vector<Real> & get_detady() const
341  { calculate_map = true; return this->_fe_map->get_detady(); }
342 
347  const std::vector<Real> & get_detadz() const
348  { calculate_map = true; return this->_fe_map->get_detadz(); }
349 
354  const std::vector<Real> & get_dzetadx() const
355  { calculate_map = true; return this->_fe_map->get_dzetadx(); }
356 
361  const std::vector<Real> & get_dzetady() const
362  { calculate_map = true; return this->_fe_map->get_dzetady(); }
363 
368  const std::vector<Real> & get_dzetadz() const
369  { calculate_map = true; return this->_fe_map->get_dzetadz(); }
370 
374  const std::vector<std::vector<Point>> & get_tangents() const
375  { calculate_map = true; return this->_fe_map->get_tangents(); }
376 
380  const std::vector<Point> & get_normals() const
381  { calculate_map = true; return this->_fe_map->get_normals(); }
382 
383 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
384 
387  const std::vector<Real> & get_curvatures() const
388  { calculate_map = true; return this->_fe_map->get_curvatures();}
389 
390 #endif
391 
396  virtual void attach_quadrature_rule (QBase * q) = 0;
397 
403  virtual unsigned int n_shape_functions () const = 0;
404 
409  virtual unsigned int n_quadrature_points () const = 0;
410 
416  ElemType get_type() const { return elem_type; }
417 
422  unsigned int get_p_level() const { return _p_level; }
423 
427  FEType get_fe_type() const { return fe_type; }
428 
432  Order get_order() const { return static_cast<Order>(fe_type.order + _p_level); }
433 
437  void set_fe_order(int new_order) { fe_type.order = new_order; }
438 
442  virtual FEContinuity get_continuity() const = 0;
443 
448  virtual bool is_hierarchic() const = 0;
449 
453  FEFamily get_family() const { return fe_type.family; }
454 
458  const FEMap & get_fe_map() const { return *_fe_map.get(); }
459  FEMap & get_fe_map() { return *_fe_map.get(); }
460 
464  void print_JxW(std::ostream & os) const;
465 
471  virtual void print_phi(std::ostream & os) const =0;
472 
478  virtual void print_dphi(std::ostream & os) const =0;
479 
480 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
481 
487  virtual void print_d2phi(std::ostream & os) const =0;
488 
489 #endif
490 
495  void print_xyz(std::ostream & os) const;
496 
500  void print_info(std::ostream & os) const;
501 
505  friend std::ostream & operator << (std::ostream & os, const FEAbstract & fe);
506 
507 
508 protected:
509 
522  virtual void compute_shape_functions(const Elem *, const std::vector<Point> & ) =0;
523 
524  std::unique_ptr<FEMap> _fe_map;
525 
526 
530  const unsigned int dim;
531 
536  mutable bool calculations_started;
537 
541  mutable bool calculate_map;
542 
546  mutable bool calculate_phi;
547 
551  mutable bool calculate_dphi;
552 
553 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
554 
557  mutable bool calculate_d2phi;
558 #else
559  // Otherwise several interfaces need to be redone.
560  const bool calculate_d2phi=false;
561 
562 #endif
563 
567  mutable bool calculate_curl_phi;
568 
572  mutable bool calculate_div_phi;
573 
577  mutable bool calculate_dphiref;
578 
579 
586 
592 
597  unsigned int _p_level;
598 
603 
609 
616  virtual bool shapes_need_reinit() const = 0;
617 
618 };
619 
620 } // namespace libMesh
621 
622 #endif // LIBMESH_FE_ABSTRACT_H
libMesh::QBase
The QBase class provides the basic functionality from which various quadrature rules can be derived.
Definition: quadrature.h:61
libMesh::FEAbstract::n_shape_functions
virtual unsigned int n_shape_functions() const =0
libMesh::FEAbstract::get_order
Order get_order() const
Definition: fe_abstract.h:432
libMesh::PeriodicBoundaries
We're using a class instead of a typedef to allow forward declarations and future flexibility.
Definition: periodic_boundaries.h:51
libMesh::FEAbstract::get_d2xyzdxideta
const std::vector< RealGradient > & get_d2xyzdxideta() const
Definition: fe_abstract.h:291
libMesh::FEAbstract::side_map
virtual void side_map(const Elem *elem, const Elem *side, const unsigned int s, const std::vector< Point > &reference_side_points, std::vector< Point > &reference_points)=0
Computes the reference space quadrature points on the side of an element based on the side quadrature...
libMesh::FEAbstract::get_xyz
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:237
libMesh::FEAbstract::calculate_dphi
bool calculate_dphi
Should we calculate shape function gradients?
Definition: fe_abstract.h:551
libMesh::FEType::family
FEFamily family
The type of finite element.
Definition: fe_type.h:203
libMesh::FEAbstract::shapes_need_reinit
virtual bool shapes_need_reinit() const =0
libMesh::FEAbstract::print_xyz
void print_xyz(std::ostream &os) const
Prints the spatial location of each quadrature point (on the physical element).
Definition: fe_abstract.C:812
libMesh::FEAbstract::get_normals
const std::vector< Point > & get_normals() const
Definition: fe_abstract.h:380
libMesh::FEAbstract::fe_type
FEType fe_type
The finite element type for this object.
Definition: fe_abstract.h:585
libMesh::ReferenceCountedObject
This class implements reference counting.
Definition: reference_counted_object.h:65
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::FEAbstract::get_tangents
const std::vector< std::vector< Point > > & get_tangents() const
Definition: fe_abstract.h:374
libMesh::FEAbstract::get_d2xyzdzeta2
const std::vector< RealGradient > & get_d2xyzdzeta2() const
Definition: fe_abstract.h:285
libMesh::FEAbstract::get_detadx
const std::vector< Real > & get_detadx() const
Definition: fe_abstract.h:333
libMesh::FEAbstract
This class forms the foundation from which generic finite elements may be derived.
Definition: fe_abstract.h:100
libMesh::FEAbstract::print_dphi
virtual void print_dphi(std::ostream &os) const =0
Prints the value of each shape function's derivative at each quadrature point.
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::FEAbstract::calculate_dphiref
bool calculate_dphiref
Should we calculate reference shape function gradients?
Definition: fe_abstract.h:577
libMesh::FEAbstract::get_JxW
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:244
libMesh::FEAbstract::calculations_started
bool calculations_started
Have calculations with this object already been started? Then all get_* functions should already have...
Definition: fe_abstract.h:536
libMesh::FEAbstract::get_dzetadx
const std::vector< Real > & get_dzetadx() const
Definition: fe_abstract.h:354
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::FEAbstract::get_dxidy
const std::vector< Real > & get_dxidy() const
Definition: fe_abstract.h:319
libMesh::FEAbstract::calculate_phi
bool calculate_phi
Should we calculate shape functions?
Definition: fe_abstract.h:546
libMesh::FEAbstract::calculate_map
bool calculate_map
Are we calculating mapping functions?
Definition: fe_abstract.h:541
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::FEAbstract::get_family
FEFamily get_family() const
Definition: fe_abstract.h:453
libMesh::FEAbstract::calculate_curl_phi
bool calculate_curl_phi
Should we calculate shape function curls?
Definition: fe_abstract.h:567
libMesh::FEAbstract::get_d2xyzdeta2
const std::vector< RealGradient > & get_d2xyzdeta2() const
Definition: fe_abstract.h:279
libMesh::FEAbstract::get_dxyzdzeta
const std::vector< RealGradient > & get_dxyzdzeta() const
Definition: fe_abstract.h:265
libMesh::FEAbstract::shapes_on_quadrature
bool shapes_on_quadrature
A flag indicating if current data structures correspond to quadrature rule points.
Definition: fe_abstract.h:608
libMesh::FEAbstract::set_fe_order
void set_fe_order(int new_order)
Sets the base FE order of the finite element.
Definition: fe_abstract.h:437
libMesh::FEAbstract::operator<<
friend std::ostream & operator<<(std::ostream &os, const FEAbstract &fe)
Same as above, but allows you to print to a stream.
Definition: fe_abstract.C:834
libMesh::FEAbstract::get_detady
const std::vector< Real > & get_detady() const
Definition: fe_abstract.h:340
libMesh::FEAbstract::get_fe_map
const FEMap & get_fe_map() const
Definition: fe_abstract.h:458
libMesh::FEAbstract::get_d2xyzdetadzeta
const std::vector< RealGradient > & get_d2xyzdetadzeta() const
Definition: fe_abstract.h:303
libMesh::FEAbstract::get_continuity
virtual FEContinuity get_continuity() const =0
libMesh::FEAbstract::print_info
void print_info(std::ostream &os) const
Prints all the relevant information about the current element.
Definition: fe_abstract.C:818
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::FEAbstract::get_p_level
unsigned int get_p_level() const
Definition: fe_abstract.h:422
libMesh::FEAbstract::get_dim
unsigned int get_dim() const
Definition: fe_abstract.h:230
libMesh::FEAbstract::print_phi
virtual void print_phi(std::ostream &os) const =0
Prints the value of each shape function at each quadrature point.
libMesh::FEAbstract::calculate_div_phi
bool calculate_div_phi
Should we calculate shape function divergences?
Definition: fe_abstract.h:572
libMesh::FEAbstract::build
static std::unique_ptr< FEAbstract > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
Definition: fe_abstract.C:72
libMesh::FEAbstract::get_dzetady
const std::vector< Real > & get_dzetady() const
Definition: fe_abstract.h:361
libMesh::FEAbstract::qrule
QBase * qrule
A pointer to the quadrature rule employed.
Definition: fe_abstract.h:602
libMesh::FEAbstract::get_type
ElemType get_type() const
Definition: fe_abstract.h:416
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::FEAbstract::_fe_map
std::unique_ptr< FEMap > _fe_map
Definition: fe_abstract.h:524
libMesh::FEAbstract::edge_reinit
virtual void edge_reinit(const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *pts=nullptr, const std::vector< Real > *weights=nullptr)=0
Reinitializes all the physical element-dependent data based on the edge of the element elem.
libMesh::FEAbstract::get_curvatures
const std::vector< Real > & get_curvatures() const
Definition: fe_abstract.h:387
libMesh::FEAbstract::get_dxidz
const std::vector< Real > & get_dxidz() const
Definition: fe_abstract.h:326
libMesh::FEAbstract::get_fe_map
FEMap & get_fe_map()
Definition: fe_abstract.h:459
libMesh::FEAbstract::n_quadrature_points
virtual unsigned int n_quadrature_points() const =0
libMesh::FEAbstract::get_dzetadz
const std::vector< Real > & get_dzetadz() const
Definition: fe_abstract.h:368
libMesh::FEAbstract::dim
const unsigned int dim
The dimensionality of the object.
Definition: fe_abstract.h:530
libMesh::FEAbstract::get_dxidx
const std::vector< Real > & get_dxidx() const
Definition: fe_abstract.h:312
libMesh::FEAbstract::compute_node_constraints
static void compute_node_constraints(NodeConstraints &constraints, const Elem *elem)
Computes the nodal constraint contributions (for non-conforming adapted meshes), using Lagrange geome...
Definition: fe_abstract.C:845
libMesh::FEAbstract::get_fe_type
FEType get_fe_type() const
Definition: fe_abstract.h:427
libMesh::FEAbstract::get_detadz
const std::vector< Real > & get_detadz() const
Definition: fe_abstract.h:347
libMesh::FEAbstract::compute_periodic_node_constraints
static void compute_periodic_node_constraints(NodeConstraints &constraints, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const Elem *elem)
Computes the node position constraint equation contributions (for meshes with periodic boundary condi...
Definition: fe_abstract.C:990
libMesh::FEAbstract::get_dxyzdeta
const std::vector< RealGradient > & get_dxyzdeta() const
Definition: fe_abstract.h:258
libMesh::FEAbstract::on_reference_element
static bool on_reference_element(const Point &p, const ElemType t, const Real eps=TOLERANCE)
Definition: fe_abstract.C:606
libMesh::NodeConstraints
The Node constraint storage format.
Definition: dof_map.h:153
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::FEType::order
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:197
libMesh::FEFamily
FEFamily
Definition: enum_fe_family.h:34
libMesh::FEAbstract::~FEAbstract
virtual ~FEAbstract()
Destructor.
Definition: fe_abstract.C:67
libMesh::FEAbstract::_p_level
unsigned int _p_level
The p refinement level the current data structures are set up for.
Definition: fe_abstract.h:597
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::FEAbstract::print_d2phi
virtual void print_d2phi(std::ostream &os) const =0
Prints the value of each shape function's second derivatives at each quadrature point.
libMesh::FEAbstract::elem_type
ElemType elem_type
The element type the current data structures are set up for.
Definition: fe_abstract.h:591
libMesh::FEAbstract::get_refspace_nodes
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:308
libMesh::FEAbstract::compute_shape_functions
virtual void compute_shape_functions(const Elem *, const std::vector< Point > &)=0
After having updated the jacobian and the transformation from local to global coordinates in FEMap::c...
libMesh::FEMap
Class contained in FE that encapsulates mapping (i.e.
Definition: fe_map.h:47
libMesh::FEAbstract::get_dxyzdxi
const std::vector< RealGradient > & get_dxyzdxi() const
Definition: fe_abstract.h:251
libMesh::FEAbstract::get_d2xyzdxi2
const std::vector< RealGradient > & get_d2xyzdxi2() const
Definition: fe_abstract.h:273
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEAbstract::reinit
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)=0
This is at the core of this class.
libMesh::FEAbstract::attach_quadrature_rule
virtual void attach_quadrature_rule(QBase *q)=0
Provides the class with the quadrature rule.
libMesh::FEContinuity
FEContinuity
Definition: enum_fe_family.h:77
libMesh::FEAbstract::FEAbstract
FEAbstract(const unsigned int dim, const FEType &fet)
Constructor.
Definition: fe_abstract.C:46
libMesh::PointLocatorBase
This is the base class for point locators.
Definition: point_locator_base.h:62
libMesh::FEAbstract::is_hierarchic
virtual bool is_hierarchic() const =0
libMesh::FEAbstract::print_JxW
void print_JxW(std::ostream &os) const
Prints the Jacobian times the weight for each quadrature point.
Definition: fe_abstract.C:805
libMesh::FEAbstract::calculate_d2phi
bool calculate_d2phi
Should we calculate shape function hessians?
Definition: fe_abstract.h:557
int
void ErrorVector unsigned int
Definition: adjoints_ex3.C:360
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::FEAbstract::get_d2xyzdxidzeta
const std::vector< RealGradient > & get_d2xyzdxidzeta() const
Definition: fe_abstract.h:297