libMesh
fe_transformation_base.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 #ifndef LIBMESH_FE_TRANSFORMATION_BASE_H
19 #define LIBMESH_FE_TRANSFORMATION_BASE_H
20 
21 #include "libmesh/fe_base.h"
22 
23 namespace libMesh
24 {
25 
26 //Forward Declarations
27 template<typename T> class FEGenericBase;
28 template<typename T> class H1FETransformation;
29 template<typename T> class HCurlFETransformation;
30 template<typename T> class HDivFETransformation;
31 class FEType;
32 
42 template<typename OutputShape>
44 {
45 public:
46 
47  FETransformationBase() = default;
48  virtual ~FETransformationBase() = default;
49 
53  static std::unique_ptr<FETransformationBase<OutputShape>> build(const FEType & type);
54 
58  virtual void init_map_phi(const FEGenericBase<OutputShape> & fe) const = 0;
59 
63  virtual void init_map_dphi(const FEGenericBase<OutputShape> & fe) const = 0;
64 
68  virtual void init_map_d2phi(const FEGenericBase<OutputShape> & fe) const = 0;
69 
74  virtual void map_phi(const unsigned int dim,
75  const Elem * const elem,
76  const std::vector<Point> & qp,
77  const FEGenericBase<OutputShape> & fe,
78  std::vector<std::vector<OutputShape>> & phi,
79  bool add_p_level = true) const = 0;
80 
85  virtual void map_dphi(const unsigned int dim,
86  const Elem * const elem,
87  const std::vector<Point> & qp,
88  const FEGenericBase<OutputShape> & fe,
89  std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> & dphi,
90  std::vector<std::vector<OutputShape>> & dphidx,
91  std::vector<std::vector<OutputShape>> & dphidy,
92  std::vector<std::vector<OutputShape>> & dphidz) const = 0;
93 
94 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
95 
99  virtual void map_d2phi(const unsigned int dim,
100  const std::vector<Point> & qp,
101  const FEGenericBase<OutputShape> & fe,
102  std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor>> & d2phi,
103  std::vector<std::vector<OutputShape>> & d2phidx2,
104  std::vector<std::vector<OutputShape>> & d2phidxdy,
105  std::vector<std::vector<OutputShape>> & d2phidxdz,
106  std::vector<std::vector<OutputShape>> & d2phidy2,
107  std::vector<std::vector<OutputShape>> & d2phidydz,
108  std::vector<std::vector<OutputShape>> & d2phidz2) const = 0;
109 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
110 
111 
116  virtual void map_curl(const unsigned int dim,
117  const Elem * const elem,
118  const std::vector<Point> & qp,
119  const FEGenericBase<OutputShape> & fe,
120  std::vector<std::vector<OutputShape>> & curl_phi) const = 0;
121 
126  virtual void map_div(const unsigned int dim,
127  const Elem * const elem,
128  const std::vector<Point> & qp,
129  const FEGenericBase<OutputShape> & fe,
130  std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence>> & div_phi) const = 0;
131 
132 }; // class FETransformationBase
133 
134 }
135 
136 #endif // LIBMESH_FE_TRANSFORMATION_BASE_H
This class handles the computation of the shape functions in the physical domain for HDiv conforming ...
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
virtual void map_div(const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputDivergence >> &div_phi) const =0
Evaluates the shape function divergence in physical coordinates based on proper finite element transf...
virtual void init_map_d2phi(const FEGenericBase< OutputShape > &fe) const =0
Pre-requests any necessary data from FEMap.
virtual void map_d2phi(const unsigned int dim, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputTensor >> &d2phi, std::vector< std::vector< OutputShape >> &d2phidx2, std::vector< std::vector< OutputShape >> &d2phidxdy, std::vector< std::vector< OutputShape >> &d2phidxdz, std::vector< std::vector< OutputShape >> &d2phidy2, std::vector< std::vector< OutputShape >> &d2phidydz, std::vector< std::vector< OutputShape >> &d2phidz2) const =0
Evaluates shape function Hessians in physical coordinates based on proper finite element transformati...
This class handles the computation of the shape functions in the physical domain for HCurl conforming...
unsigned int dim
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
The libMesh namespace provides an interface to certain functionality in the library.
TensorTools::DecrementRank< OutputShape >::type OutputDivergence
Definition: fe_base.h:122
virtual void map_phi(const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< OutputShape >> &phi, bool add_p_level=true) const =0
Evaluates shape functions in physical coordinates based on proper finite element transformation.
This class handles the computation of the shape functions in the physical domain for H1 conforming el...
static std::unique_ptr< FETransformationBase< OutputShape > > build(const FEType &type)
Builds an FETransformation object based on the finite element type.
virtual ~FETransformationBase()=default
virtual void map_dphi(const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputGradient >> &dphi, std::vector< std::vector< OutputShape >> &dphidx, std::vector< std::vector< OutputShape >> &dphidy, std::vector< std::vector< OutputShape >> &dphidz) const =0
Evaluates shape function gradients in physical coordinates based on proper finite element transformat...
virtual void init_map_phi(const FEGenericBase< OutputShape > &fe) const =0
Pre-requests any necessary data from FEMap.
virtual void init_map_dphi(const FEGenericBase< OutputShape > &fe) const =0
Pre-requests any necessary data from FEMap.
This class handles the computation of the shape functions in the physical domain. ...
Definition: fe_base.h:54
virtual void map_curl(const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< OutputShape >> &curl_phi) const =0
Evaluates the shape function curl in physical coordinates based on proper finite element transformati...
This class forms the foundation from which generic finite elements may be derived.