libMesh
Public Member Functions | Static Public Member Functions | List of all members
libMesh::HDivFETransformation< OutputShape > Class Template Reference

This class handles the computation of the shape functions in the physical domain for HDiv conforming elements. More...

#include <fe_transformation_base.h>

Inheritance diagram for libMesh::HDivFETransformation< OutputShape >:
[legend]

Public Member Functions

 HDivFETransformation ()
 
virtual ~HDivFETransformation ()=default
 
virtual void init_map_phi (const FEGenericBase< OutputShape > &fe) const override
 Pre-requests any necessary data from FEMap. More...
 
virtual void init_map_dphi (const FEGenericBase< OutputShape > &fe) const override
 Pre-requests any necessary data from FEMap. More...
 
virtual void init_map_d2phi (const FEGenericBase< OutputShape > &fe) const override
 Pre-requests any necessary data from FEMap. More...
 
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 override
 Evaluates shape functions in physical coordinates for \( H(div) \) conforming elements. More...
 
virtual void map_dphi (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< OutputShape > &, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputGradient >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &) const override
 Evaluates shape function gradients in physical coordinates for \( H(div) \) conforming elements. More...
 
virtual void map_d2phi (const unsigned int, const std::vector< Point > &, const FEGenericBase< OutputShape > &, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputTensor >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &, std::vector< std::vector< OutputShape >> &) const override
 Evaluates shape function Hessians in physical coordinates based on \( H(div) \) conforming finite element transformation. More...
 
virtual void map_curl (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< OutputShape > &, std::vector< std::vector< OutputShape >> &) const override
 Evaluates the shape function curl in physical coordinates based on \( H(div) \) conforming finite element transformation. More...
 
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 override
 Evaluates the shape function divergence in physical coordinates based on \( H(div) \) conforming finite element transformation. More...
 
template<>
void init_map_phi (const FEGenericBase< Real > &) const
 
template<>
void init_map_dphi (const FEGenericBase< Real > &) const
 
template<>
void init_map_d2phi (const FEGenericBase< Real > &) const
 
template<>
void map_phi (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< Real > &, std::vector< std::vector< Real >> &, bool) const
 
template<>
void map_div (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< Real > &, std::vector< std::vector< Real >> &) const
 

Static Public Member Functions

static std::unique_ptr< FETransformationBase< OutputShape > > build (const FEType &type)
 Builds an FETransformation object based on the finite element type. More...
 

Detailed Description

template<typename OutputShape>
class libMesh::HDivFETransformation< OutputShape >

This class handles the computation of the shape functions in the physical domain for HDiv conforming elements.

This class assumes the FEGenericBase object has been initialized in the reference domain (i.e. init_shape_functions has been called).

Author
Nuno Nobre
Date
2023

Definition at line 30 of file fe_transformation_base.h.

Constructor & Destructor Documentation

◆ HDivFETransformation()

template<typename OutputShape >
libMesh::HDivFETransformation< OutputShape >::HDivFETransformation ( )
inline

Definition at line 40 of file hdiv_fe_transformation.h.

41  : FETransformationBase<OutputShape>() {}

◆ ~HDivFETransformation()

template<typename OutputShape >
virtual libMesh::HDivFETransformation< OutputShape >::~HDivFETransformation ( )
virtualdefault

Member Function Documentation

◆ build()

template<typename OutputShape >
std::unique_ptr< FETransformationBase< OutputShape > > libMesh::FETransformationBase< OutputShape >::build ( const FEType type)
staticinherited

Builds an FETransformation object based on the finite element type.

Definition at line 32 of file fe_transformation_base.C.

References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::Utility::enum_to_string(), libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::HIERARCHIC_VEC, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::L2_HIERARCHIC, libMesh::L2_HIERARCHIC_VEC, libMesh::L2_LAGRANGE, libMesh::L2_LAGRANGE_VEC, libMesh::L2_RAVIART_THOMAS, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::MONOMIAL_VEC, libMesh::NEDELEC_ONE, libMesh::RATIONAL_BERNSTEIN, libMesh::RAVIART_THOMAS, libMesh::SCALAR, libMesh::SIDE_HIERARCHIC, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.

33 {
34  switch (fe_type.family)
35  {
36  // H1 Conforming Elements
37  case LAGRANGE:
38  case HIERARCHIC:
39  case BERNSTEIN:
40  case SZABAB:
41  case CLOUGH: // PB: Really H2
42  case HERMITE: // PB: Really H2
43  case SUBDIVISION:
44  case LAGRANGE_VEC:
45  case HIERARCHIC_VEC:
46  case MONOMIAL: // PB: Shouldn't this be L2 conforming?
47  case MONOMIAL_VEC: // PB: Shouldn't this be L2 conforming?
48  case XYZ: // PB: Shouldn't this be L2 conforming?
49  case RATIONAL_BERNSTEIN:
50  case L2_HIERARCHIC:
51  case L2_HIERARCHIC_VEC:
52  case SIDE_HIERARCHIC:
53  case L2_LAGRANGE: // PB: Shouldn't this be L2 conforming?
54  case L2_LAGRANGE_VEC: // PB: Shouldn't this be L2 conforming?
55  case JACOBI_20_00: // PB: For infinite elements...
56  case JACOBI_30_00: // PB: For infinite elements...
57  return std::make_unique<H1FETransformation<OutputShape>>();
58 
59  // HCurl Conforming Elements
60  case NEDELEC_ONE:
61  return std::make_unique<HCurlFETransformation<OutputShape>>();
62 
63  // HDiv Conforming Elements
64  case RAVIART_THOMAS:
65  case L2_RAVIART_THOMAS:
66  return std::make_unique<HDivFETransformation<OutputShape>>();
67 
68  // L2 Conforming Elements
69 
70  // Other...
71  case SCALAR:
72  // Should never need this for SCALARs
73  return std::make_unique<H1FETransformation<OutputShape>>();
74 
75  default:
76  libmesh_error_msg("Unknown family = " << Utility::enum_to_string(fe_type.family));
77  }
78 }
std::string enum_to_string(const T e)

◆ init_map_d2phi() [1/2]

template<typename OutputShape >
void libMesh::HDivFETransformation< OutputShape >::init_map_d2phi ( const FEGenericBase< OutputShape > &  fe) const
overridevirtual

Pre-requests any necessary data from FEMap.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 42 of file hdiv_fe_transformation.C.

References libMesh::FEMap::get_d2xidxyz2(), libMesh::FEMap::get_dxidx(), and libMesh::FEAbstract::get_fe_map().

43 {
44  fe.get_fe_map().get_dxidx();
45 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
46  fe.get_fe_map().get_d2xidxyz2();
47 #endif
48 }

◆ init_map_d2phi() [2/2]

template<>
void libMesh::HDivFETransformation< Real >::init_map_d2phi ( const FEGenericBase< Real > &  ) const

Definition at line 209 of file hdiv_fe_transformation.C.

210 {
211  libmesh_error_msg("HDiv transformations only make sense for vector-valued elements.");
212 }

◆ init_map_dphi() [1/2]

template<typename OutputShape >
void libMesh::HDivFETransformation< OutputShape >::init_map_dphi ( const FEGenericBase< OutputShape > &  fe) const
overridevirtual

Pre-requests any necessary data from FEMap.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 34 of file hdiv_fe_transformation.C.

References libMesh::FEMap::get_dxidx(), and libMesh::FEAbstract::get_fe_map().

35 {
36  fe.get_fe_map().get_dxidx();
37 }

◆ init_map_dphi() [2/2]

template<>
void libMesh::HDivFETransformation< Real >::init_map_dphi ( const FEGenericBase< Real > &  ) const

Definition at line 203 of file hdiv_fe_transformation.C.

204 {
205  libmesh_error_msg("HDiv transformations only make sense for vector-valued elements.");
206 }

◆ init_map_phi() [1/2]

template<typename OutputShape >
void libMesh::HDivFETransformation< OutputShape >::init_map_phi ( const FEGenericBase< OutputShape > &  fe) const
overridevirtual

Pre-requests any necessary data from FEMap.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 26 of file hdiv_fe_transformation.C.

References libMesh::FEMap::get_dxidx(), and libMesh::FEAbstract::get_fe_map().

27 {
28  fe.get_fe_map().get_dxidx();
29 }

◆ init_map_phi() [2/2]

template<>
void libMesh::HDivFETransformation< Real >::init_map_phi ( const FEGenericBase< Real > &  ) const

Definition at line 197 of file hdiv_fe_transformation.C.

198 {
199  libmesh_error_msg("HDiv transformations only make sense for vector-valued elements.");
200 }

◆ map_curl()

template<typename OutputShape >
virtual void libMesh::HDivFETransformation< OutputShape >::map_curl ( const unsigned int  ,
const Elem const,
const std::vector< Point > &  ,
const FEGenericBase< OutputShape > &  ,
std::vector< std::vector< OutputShape >> &   
) const
inlineoverridevirtual

Evaluates the shape function curl in physical coordinates based on \( H(div) \) conforming finite element transformation.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 115 of file hdiv_fe_transformation.h.

120  {
121  libmesh_warning("WARNING: Shape function curls for HDiv elements are not currently being computed!");
122  }

◆ map_d2phi()

template<typename OutputShape >
virtual void libMesh::HDivFETransformation< OutputShape >::map_d2phi ( const unsigned int  ,
const std::vector< Point > &  ,
const FEGenericBase< OutputShape > &  ,
std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputTensor >> &  ,
std::vector< std::vector< OutputShape >> &  ,
std::vector< std::vector< OutputShape >> &  ,
std::vector< std::vector< OutputShape >> &  ,
std::vector< std::vector< OutputShape >> &  ,
std::vector< std::vector< OutputShape >> &  ,
std::vector< std::vector< OutputShape >> &   
) const
inlineoverridevirtual

Evaluates shape function Hessians in physical coordinates based on \( H(div) \) conforming finite element transformation.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 96 of file hdiv_fe_transformation.h.

106  {
107  libmesh_warning("WARNING: Shape function Hessians for HDiv elements are not currently being computed!");
108  }

◆ map_div() [1/2]

template<typename OutputShape >
void libMesh::HDivFETransformation< OutputShape >::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
overridevirtual

Evaluates the shape function divergence in physical coordinates based on \( H(div) \) conforming finite element transformation.

The transformation is \( \nabla \cdot \phi = J^{-1} \nabla \cdot \hat{\phi} \) where \( J = \det( dx/d\xi ) \)

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 143 of file hdiv_fe_transformation.C.

References dim, libMesh::FEGenericBase< OutputType >::get_dphideta(), libMesh::FEGenericBase< OutputType >::get_dphidxi(), libMesh::FEGenericBase< OutputType >::get_dphidzeta(), libMesh::FEAbstract::get_fe_map(), libMesh::FEMap::get_jacobian(), and libMesh::index_range().

148 {
149  switch (dim)
150  {
151  case 0:
152  case 1:
153  libmesh_error_msg("These element transformations only make sense in 2D and 3D.");
154 
155  case 2:
156  {
157  const std::vector<std::vector<OutputShape>> & dphi_dxi = fe.get_dphidxi();
158  const std::vector<std::vector<OutputShape>> & dphi_deta = fe.get_dphideta();
159 
160  const std::vector<Real> & J = fe.get_fe_map().get_jacobian();
161 
162  // div(phi) = J^{-1} * div(\hat{phi})
163  for (auto i : index_range(div_phi))
164  for (auto p : index_range(div_phi[i]))
165  {
166  div_phi[i][p] = (dphi_dxi[i][p](0) + dphi_deta[i][p](1))/J[p];
167  }
168 
169  break;
170  }
171  case 3:
172  {
173  const std::vector<std::vector<OutputShape>> & dphi_dxi = fe.get_dphidxi();
174  const std::vector<std::vector<OutputShape>> & dphi_deta = fe.get_dphideta();
175  const std::vector<std::vector<OutputShape>> & dphi_dzeta = fe.get_dphidzeta();
176 
177  const std::vector<Real> & J = fe.get_fe_map().get_jacobian();
178 
179  // div(phi) = J^{-1} * div(\hat{phi})
180  for (auto i : index_range(div_phi))
181  for (auto p : index_range(div_phi[i]))
182  {
183  div_phi[i][p] = (dphi_dxi[i][p](0) + dphi_deta[i][p](1) + dphi_dzeta[i][p](2))/J[p];
184  }
185 
186  break;
187  }
188 
189  default:
190  libmesh_error_msg("Invalid dim = " << dim);
191  } // switch(dim)
192 }
unsigned int dim
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117

◆ map_div() [2/2]

template<>
void libMesh::HDivFETransformation< Real >::map_div ( const unsigned int  ,
const Elem const,
const std::vector< Point > &  ,
const FEGenericBase< Real > &  ,
std::vector< std::vector< Real >> &   
) const

Definition at line 226 of file hdiv_fe_transformation.C.

231 {
232  libmesh_error_msg("HDiv transformations only make sense for vector-valued elements.");
233 }

◆ map_dphi()

template<typename OutputShape >
virtual void libMesh::HDivFETransformation< OutputShape >::map_dphi ( const unsigned int  ,
const Elem const,
const std::vector< Point > &  ,
const FEGenericBase< OutputShape > &  ,
std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputGradient >> &  ,
std::vector< std::vector< OutputShape >> &  ,
std::vector< std::vector< OutputShape >> &  ,
std::vector< std::vector< OutputShape >> &   
) const
inlineoverridevirtual

Evaluates shape function gradients in physical coordinates for \( H(div) \) conforming elements.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 79 of file hdiv_fe_transformation.h.

87  {
88  libmesh_warning("WARNING: Shape function gradients for HDiv elements are not currently being computed!");
89  }

◆ map_phi() [1/2]

template<typename OutputShape >
void libMesh::HDivFETransformation< OutputShape >::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
overridevirtual

Evaluates shape functions in physical coordinates for \( H(div) \) conforming elements.

In this case \( \phi = J^{-1} (dx/d\xi) \hat{\phi} \), where \( (dx/d\xi) \) is the Jacobian matrix of the element map and J = ( dx/d ).

Note
Here \( x, \xi \) are vectors.

Implements libMesh::FETransformationBase< OutputShape >.

Definition at line 53 of file hdiv_fe_transformation.C.

References dim, libMesh::FEMap::get_dxyzdeta(), libMesh::FEMap::get_dxyzdxi(), libMesh::FEMap::get_dxyzdzeta(), libMesh::FEAbstract::get_fe_map(), libMesh::FEAbstract::get_fe_type(), libMesh::FEMap::get_jacobian(), libMesh::index_range(), libMesh::Real, and libMesh::FEInterface::shape().

59 {
60  switch (dim)
61  {
62  case 0:
63  case 1:
64  libmesh_error_msg("These element transformations only make sense in 2D and 3D.");
65 
66  case 2:
67  {
68  const std::vector<RealGradient> & dxyz_dxi = fe.get_fe_map().get_dxyzdxi();
69  const std::vector<RealGradient> & dxyz_deta = fe.get_fe_map().get_dxyzdeta();
70 
71  const std::vector<Real> & J = fe.get_fe_map().get_jacobian();
72 
73  // phi = J^{-1} * (dx/dxi) * \hat{phi}
74  for (auto i : index_range(phi))
75  for (auto p : index_range(phi[i]))
76  {
77  Real dx_dxi = dxyz_dxi[p](0);
78  Real dx_deta = dxyz_deta[p](0);
79 
80  Real dy_dxi = dxyz_dxi[p](1);
81  Real dy_deta = dxyz_deta[p](1);
82 
83  Real dz_dxi = dxyz_dxi[p](2);
84  Real dz_deta = dxyz_deta[p](2);
85 
86  // Need to temporarily cache reference shape functions
87  // We are computing mapping basis functions, so we explicitly ignore
88  // any non-zero p_level() the Elem might have.
89  OutputShape phi_ref;
90  FEInterface::shape(fe.get_fe_type(), /*extra_order=*/0, elem, i, qp[p], phi_ref);
91 
92  phi[i][p](0) = (dx_dxi*phi_ref(0) + dx_deta*phi_ref(1))/J[p];
93  phi[i][p](1) = (dy_dxi*phi_ref(0) + dy_deta*phi_ref(1))/J[p];
94  phi[i][p](2) = (dz_dxi*phi_ref(0) + dz_deta*phi_ref(1))/J[p];
95  }
96 
97  break;
98  }
99  case 3:
100  {
101  const std::vector<RealGradient> & dxyz_dxi = fe.get_fe_map().get_dxyzdxi();
102  const std::vector<RealGradient> & dxyz_deta = fe.get_fe_map().get_dxyzdeta();
103  const std::vector<RealGradient> & dxyz_dzeta = fe.get_fe_map().get_dxyzdzeta();
104 
105  const std::vector<Real> & J = fe.get_fe_map().get_jacobian();
106 
107  // phi = J^{-1} * (dx/dxi) * \hat{phi}
108  for (auto i : index_range(phi))
109  for (auto p : index_range(phi[i]))
110  {
111  Real dx_dxi = dxyz_dxi[p](0);
112  Real dx_deta = dxyz_deta[p](0);
113  Real dx_dzeta = dxyz_dzeta[p](0);
114 
115  Real dy_dxi = dxyz_dxi[p](1);
116  Real dy_deta = dxyz_deta[p](1);
117  Real dy_dzeta = dxyz_dzeta[p](1);
118 
119  Real dz_dxi = dxyz_dxi[p](2);
120  Real dz_deta = dxyz_deta[p](2);
121  Real dz_dzeta = dxyz_dzeta[p](2);
122 
123  // Need to temporarily cache reference shape functions
124  // We are computing mapping basis functions, so we explicitly ignore
125  // any non-zero p_level() the Elem might have.
126  OutputShape phi_ref;
127  FEInterface::shape(fe.get_fe_type(), /*extra_order=*/0, elem, i, qp[p], phi_ref);
128 
129  phi[i][p](0) = (dx_dxi*phi_ref(0) + dx_deta*phi_ref(1) + dx_dzeta*phi_ref(2))/J[p];
130  phi[i][p](1) = (dy_dxi*phi_ref(0) + dy_deta*phi_ref(1) + dy_dzeta*phi_ref(2))/J[p];
131  phi[i][p](2) = (dz_dxi*phi_ref(0) + dz_deta*phi_ref(1) + dz_dzeta*phi_ref(2))/J[p];
132  }
133 
134  break;
135  }
136 
137  default:
138  libmesh_error_msg("Invalid dim = " << dim);
139  } // switch(dim)
140 }
unsigned int dim
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
Definition: fe_interface.C:760
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117

◆ map_phi() [2/2]

template<>
void libMesh::HDivFETransformation< Real >::map_phi ( const unsigned int  ,
const Elem const,
const std::vector< Point > &  ,
const FEGenericBase< Real > &  ,
std::vector< std::vector< Real >> &  ,
bool   
) const

Definition at line 215 of file hdiv_fe_transformation.C.

221 {
222  libmesh_error_msg("HDiv transformations only make sense for vector-valued elements.");
223 }

The documentation for this class was generated from the following files: