1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 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 class FEType;
31 
32 /**
33  * This class handles the computation of the shape functions in the physical domain.
34  * Derived classes implement the particular mapping: H1, HCurl, HDiv, or L2. This
35  * class assumes the \p FEGenericBase object has been initialized in the reference
36  * domain (i.e. init_shape_functions has been called).
37  *
38  * \author Paul T. Bauman
39  * \date 2012
40  */
41 template<typename OutputShape>
42 class FETransformationBase
43 {
44 public:
45 
FETransformationBase()46   FETransformationBase(){}
~FETransformationBase()47   virtual ~FETransformationBase(){}
48 
49   /**
50    * Builds an FETransformation object based on the finite element type
51    */
52   static std::unique_ptr<FETransformationBase<OutputShape>> build(const FEType & type);
53 
54   /**
55    * Pre-requests any necessary data from FEMap
56    */
57   virtual void init_map_phi(const FEGenericBase<OutputShape> & fe) const = 0;
58 
59   /**
60    * Pre-requests any necessary data from FEMap
61    */
62   virtual void init_map_dphi(const FEGenericBase<OutputShape> & fe) const = 0;
63 
64   /**
65    * Pre-requests any necessary data from FEMap
66    */
67   virtual void init_map_d2phi(const FEGenericBase<OutputShape> & fe) const = 0;
68 
69   /**
70    * Evaluates shape functions in physical coordinates based on proper
71    * finite element transformation.
72    */
73   virtual void map_phi(const unsigned int dim,
74                        const Elem * const elem,
75                        const std::vector<Point> & qp,
76                        const FEGenericBase<OutputShape> & fe,
77                        std::vector<std::vector<OutputShape>> & phi) const = 0;
78 
79   /**
80    * Evaluates shape function gradients in physical coordinates based on proper
81    * finite element transformation.
82    */
83   virtual void map_dphi(const unsigned int dim,
84                         const Elem * const elem,
85                         const std::vector<Point> & qp,
86                         const FEGenericBase<OutputShape> & fe,
87                         std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> & dphi,
88                         std::vector<std::vector<OutputShape>> & dphidx,
89                         std::vector<std::vector<OutputShape>> & dphidy,
90                         std::vector<std::vector<OutputShape>> & dphidz) const = 0;
91 
92 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
93   /**
94    * Evaluates shape function Hessians in physical coordinates based on proper
95    * finite element transformation.
96    */
97   virtual void map_d2phi(const unsigned int dim,
98                          const std::vector<Point> & qp,
99                          const FEGenericBase<OutputShape> & fe,
100                          std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor>> & d2phi,
101                          std::vector<std::vector<OutputShape>> & d2phidx2,
102                          std::vector<std::vector<OutputShape>> & d2phidxdy,
103                          std::vector<std::vector<OutputShape>> & d2phidxdz,
104                          std::vector<std::vector<OutputShape>> & d2phidy2,
105                          std::vector<std::vector<OutputShape>> & d2phidydz,
106                          std::vector<std::vector<OutputShape>> & d2phidz2) const = 0;
107 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
108 
109 
110   /**
111    * Evaluates the shape function curl in physical coordinates based on proper
112    * finite element transformation.
113    */
114   virtual void map_curl(const unsigned int dim,
115                         const Elem * const elem,
116                         const std::vector<Point> & qp,
117                         const FEGenericBase<OutputShape> & fe,
118                         std::vector<std::vector<OutputShape>> & curl_phi) const = 0;
119 
120   /**
121    * Evaluates the shape function divergence in physical coordinates based on proper
122    * finite element transformation.
123    */
124   virtual void map_div(const unsigned int dim,
125                        const Elem * const elem,
126                        const std::vector<Point> & qp,
127                        const FEGenericBase<OutputShape> & fe,
128                        std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence>> & div_phi) const = 0;
129 
130 }; // class FETransformationBase
131 
132 }
133 
134 #endif // LIBMESH_FE_TRANSFORMATION_BASE_H
135