libMesh
hdiv_fe_transformation.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_HDIV_FE_TRANSFORMATION_H
19 #define LIBMESH_HDIV_FE_TRANSFORMATION_H
20 
21 #include "libmesh/fe_transformation_base.h"
22 
23 namespace libMesh
24 {
25 
35 template<typename OutputShape>
36 class HDivFETransformation : public FETransformationBase<OutputShape>
37 {
38 public:
39 
41  : FETransformationBase<OutputShape>() {}
42 
43  virtual ~HDivFETransformation() = default;
44 
48  virtual void init_map_phi(const FEGenericBase<OutputShape> & fe) const override;
49 
53  virtual void init_map_dphi(const FEGenericBase<OutputShape> & fe) const override;
54 
58  virtual void init_map_d2phi(const FEGenericBase<OutputShape> & fe) const override;
59 
68  virtual void map_phi(const unsigned int dim,
69  const Elem * const elem,
70  const std::vector<Point> & qp,
71  const FEGenericBase<OutputShape> & fe,
72  std::vector<std::vector<OutputShape>> & phi,
73  bool add_p_level = true) const override;
74 
79  virtual void map_dphi(const unsigned int /*dim*/,
80  const Elem * const /*elem*/,
81  const std::vector<Point> & /*qp*/,
82  const FEGenericBase<OutputShape> & /*fe*/,
83  std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> & /*dphi*/,
84  std::vector<std::vector<OutputShape>> & /*dphidx*/,
85  std::vector<std::vector<OutputShape>> & /*dphidy*/,
86  std::vector<std::vector<OutputShape>> & /*dphidz*/) const override
87  {
88  libmesh_warning("WARNING: Shape function gradients for HDiv elements are not currently being computed!");
89  }
90 
91 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
92 
96  virtual void map_d2phi(const unsigned int /*dim*/,
97  const std::vector<Point> & /*qp*/,
98  const FEGenericBase<OutputShape> & /*fe*/,
99  std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor>> & /*d2phi*/,
100  std::vector<std::vector<OutputShape>> & /*d2phidx2*/,
101  std::vector<std::vector<OutputShape>> & /*d2phidxdy*/,
102  std::vector<std::vector<OutputShape>> & /*d2phidxdz*/,
103  std::vector<std::vector<OutputShape>> & /*d2phidy2*/,
104  std::vector<std::vector<OutputShape>> & /*d2phidydz*/,
105  std::vector<std::vector<OutputShape>> & /*d2phidz2*/) const override
106  {
107  libmesh_warning("WARNING: Shape function Hessians for HDiv elements are not currently being computed!");
108  }
109 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
110 
115  virtual void map_curl(const unsigned int /*dim*/,
116  const Elem * const /*elem*/,
117  const std::vector<Point> & /*qp*/,
118  const FEGenericBase<OutputShape> & /*fe*/,
119  std::vector<std::vector<OutputShape>> & /*curl_phi*/) const override
120  {
121  libmesh_warning("WARNING: Shape function curls for HDiv elements are not currently being computed!");
122  }
123 
130  virtual void map_div(const unsigned int dim,
131  const Elem * const elem,
132  const std::vector<Point> & qp,
133  const FEGenericBase<OutputShape> & fe,
134  std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence>> & div_phi) const override;
135 
136 }; // class HDivFETransformation
137 
138 }
139 
140 #endif // LIBMESH_HDIV_FE_TRANSFORMATION_H
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 conforming finite element transfo...
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 conforming finite element t...
unsigned int dim
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
virtual void init_map_phi(const FEGenericBase< OutputShape > &fe) const override
Pre-requests any necessary data from FEMap.
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 override
Evaluates shape functions in physical coordinates for conforming elements.
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 conforming finite element transfo...
virtual void init_map_dphi(const FEGenericBase< OutputShape > &fe) const override
Pre-requests any necessary data from FEMap.
virtual void init_map_d2phi(const FEGenericBase< OutputShape > &fe) const override
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 ~HDivFETransformation()=default
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 conforming elements.
This class forms the foundation from which generic finite elements may be derived.