libMesh
hcurl_fe_transformation.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 #include "libmesh/hcurl_fe_transformation.h"
19 #include "libmesh/fe_interface.h"
20 #include "libmesh/int_range.h"
21 
22 namespace libMesh
23 {
24 
25 template<typename OutputShape>
27 {
28  fe.get_fe_map().get_dxidx();
29 }
30 
31 
32 
33 template<typename OutputShape>
35 {
36  fe.get_fe_map().get_dxidx();
37 }
38 
39 
40 
41 template<typename OutputShape>
43 {
44  fe.get_fe_map().get_dxidx();
45 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
47 #endif
48 }
49 
50 
51 
52 template<typename OutputShape>
54  const Elem * const elem,
55  const std::vector<Point> & qp,
56  const FEGenericBase<OutputShape> & fe,
57  std::vector<std::vector<OutputShape>> & phi,
58  const bool /*add_p_level*/) const
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<Real> & dxidx_map = fe.get_fe_map().get_dxidx();
69  const std::vector<Real> & dxidy_map = fe.get_fe_map().get_dxidy();
70 
71  const std::vector<Real> & detadx_map = fe.get_fe_map().get_detadx();
72  const std::vector<Real> & detady_map = fe.get_fe_map().get_detady();
73 
74  // FIXME: Need to update for 2D elements in 3D space
75  // phi = (dx/dxi)^-T * \hat{phi}
76  // In 2D:
77  // (dx/dxi)^{-1} = [ dxi/dx dxi/dy
78  // deta/dx deta/dy ]
79  //
80  // so: dxi/dx^{-T} * \hat{phi} = [ dxi/dx deta/dx [ \hat{phi}_xi
81  // dxi/dy deta/dy ] \hat{phi}_eta ]
82  //
83  // or in indicial notation: phi_j = xi_{i,j}*\hat{phi}_i
84  for (auto i : index_range(phi))
85  for (auto p : index_range(phi[i]))
86  {
87  // Need to temporarily cache reference shape functions
88  // TODO: PB: Might be worth trying to build phi_ref separately to see
89  // if we can get vectorization
90  // We are computing mapping basis functions, so we explicitly ignore
91  // any non-zero p_level() the Elem might have.
92  OutputShape phi_ref;
93  FEInterface::shape(fe.get_fe_type(), /*extra_order=*/0, elem, i, qp[p], phi_ref);
94 
95  phi[i][p](0) = dxidx_map[p]*phi_ref.slice(0) + detadx_map[p]*phi_ref.slice(1);
96 
97  phi[i][p](1) = dxidy_map[p]*phi_ref.slice(0) + detady_map[p]*phi_ref.slice(1);
98  }
99 
100  break;
101  }
102 
103  case 3:
104  {
105  const std::vector<Real> & dxidx_map = fe.get_fe_map().get_dxidx();
106  const std::vector<Real> & dxidy_map = fe.get_fe_map().get_dxidy();
107  const std::vector<Real> & dxidz_map = fe.get_fe_map().get_dxidz();
108 
109  const std::vector<Real> & detadx_map = fe.get_fe_map().get_detadx();
110  const std::vector<Real> & detady_map = fe.get_fe_map().get_detady();
111  const std::vector<Real> & detadz_map = fe.get_fe_map().get_detadz();
112 
113  const std::vector<Real> & dzetadx_map = fe.get_fe_map().get_dzetadx();
114  const std::vector<Real> & dzetady_map = fe.get_fe_map().get_dzetady();
115  const std::vector<Real> & dzetadz_map = fe.get_fe_map().get_dzetadz();
116 
117  // phi = (dx/dxi)^-T * \hat{phi}
118  // In 3D:
119  // dx/dxi^-1 = [ dxi/dx dxi/dy dxi/dz
120  // deta/dx deta/dy deta/dz
121  // dzeta/dx dzeta/dy dzeta/dz]
122  //
123  // so: dxi/dx^-T * \hat{phi} = [ dxi/dx deta/dx dzeta/dx [ \hat{phi}_xi
124  // dxi/dy deta/dy dzeta/dy \hat{phi}_eta
125  // dxi/dz deta/dz dzeta/dz ] \hat{phi}_zeta ]
126  //
127  // or in indicial notation: phi_j = xi_{i,j}*\hat{phi}_i
128  for (auto i : index_range(phi))
129  for (auto p : index_range(phi[i]))
130  {
131  // Need to temporarily cache reference shape functions
132  // TODO: PB: Might be worth trying to build phi_ref separately to see
133  // if we can get vectorization
134  // We are computing mapping basis functions, so we explicitly ignore
135  // any non-zero p_level() the Elem might have.
136  OutputShape phi_ref;
137  FEInterface::shape(fe.get_fe_type(), /*extra_order=*/0, elem, i, qp[p], phi_ref);
138 
139  phi[i][p].slice(0) = dxidx_map[p]*phi_ref.slice(0) + detadx_map[p]*phi_ref.slice(1)
140  + dzetadx_map[p]*phi_ref.slice(2);
141 
142  phi[i][p].slice(1) = dxidy_map[p]*phi_ref.slice(0) + detady_map[p]*phi_ref.slice(1)
143  + dzetady_map[p]*phi_ref.slice(2);
144 
145  phi[i][p].slice(2) = dxidz_map[p]*phi_ref.slice(0) + detadz_map[p]*phi_ref.slice(1)
146  + dzetadz_map[p]*phi_ref.slice(2);
147  }
148  break;
149  }
150 
151  default:
152  libmesh_error_msg("Invalid dim = " << dim);
153  } // switch(dim)
154 }
155 
156 template<typename OutputShape>
158  const Elem * const,
159  const std::vector<Point> &,
160  const FEGenericBase<OutputShape> & fe,
161  std::vector<std::vector<OutputShape>> & curl_phi) const
162 {
163  switch (dim)
164  {
165  case 0:
166  case 1:
167  libmesh_error_msg("These element transformations only make sense in 2D and 3D.");
168 
169  case 2:
170  {
171  const std::vector<std::vector<OutputShape>> & dphi_dxi = fe.get_dphidxi();
172  const std::vector<std::vector<OutputShape>> & dphi_deta = fe.get_dphideta();
173 
174  const std::vector<Real> & J = fe.get_fe_map().get_jacobian();
175 
176  // FIXME: I don't think this is valid for 2D elements in 3D space
177  // In 2D: curl(phi) = J^{-1} * curl(\hat{phi})
178  for (auto i : index_range(curl_phi))
179  for (auto p : index_range(curl_phi[i]))
180  {
181  curl_phi[i][p].slice(0) = curl_phi[i][p].slice(1) = 0.0;
182  curl_phi[i][p].slice(2) = ( dphi_dxi[i][p].slice(1) - dphi_deta[i][p].slice(0) )/J[p];
183  }
184 
185  break;
186  }
187 
188  case 3:
189  {
190  const std::vector<std::vector<OutputShape>> & dphi_dxi = fe.get_dphidxi();
191  const std::vector<std::vector<OutputShape>> & dphi_deta = fe.get_dphideta();
192  const std::vector<std::vector<OutputShape>> & dphi_dzeta = fe.get_dphidzeta();
193 
194  const std::vector<RealGradient> & dxyz_dxi = fe.get_fe_map().get_dxyzdxi();
195  const std::vector<RealGradient> & dxyz_deta = fe.get_fe_map().get_dxyzdeta();
196  const std::vector<RealGradient> & dxyz_dzeta = fe.get_fe_map().get_dxyzdzeta();
197 
198  const std::vector<Real> & J = fe.get_fe_map().get_jacobian();
199 
200  for (auto i : index_range(curl_phi))
201  for (auto p : index_range(curl_phi[i]))
202  {
203  Real dx_dxi = dxyz_dxi[p](0);
204  Real dx_deta = dxyz_deta[p](0);
205  Real dx_dzeta = dxyz_dzeta[p](0);
206 
207  Real dy_dxi = dxyz_dxi[p](1);
208  Real dy_deta = dxyz_deta[p](1);
209  Real dy_dzeta = dxyz_dzeta[p](1);
210 
211  Real dz_dxi = dxyz_dxi[p](2);
212  Real dz_deta = dxyz_deta[p](2);
213  Real dz_dzeta = dxyz_dzeta[p](2);
214 
215  const Real inv_jac = 1.0/J[p];
216 
217  /* In 3D: curl(phi) = J^{-1} dx/dxi * curl(\hat{phi})
218 
219  dx/dxi = [ dx/dxi dx/deta dx/dzeta
220  dy/dxi dy/deta dy/dzeta
221  dz/dxi dz/deta dz/dzeta ]
222 
223  curl(u) = [ du_z/deta - du_y/dzeta
224  du_x/dzeta - du_z/dxi
225  du_y/dxi - du_x/deta ]
226  */
227  curl_phi[i][p].slice(0) = inv_jac*( dx_dxi*( dphi_deta[i][p].slice(2) -
228  dphi_dzeta[i][p].slice(1) ) +
229  dx_deta*( dphi_dzeta[i][p].slice(0) -
230  dphi_dxi[i][p].slice(2) ) +
231  dx_dzeta*( dphi_dxi[i][p].slice(1) -
232  dphi_deta[i][p].slice(0) ) );
233 
234  curl_phi[i][p].slice(1) = inv_jac*( dy_dxi*( dphi_deta[i][p].slice(2) -
235  dphi_dzeta[i][p].slice(1) ) +
236  dy_deta*( dphi_dzeta[i][p].slice(0)-
237  dphi_dxi[i][p].slice(2) ) +
238  dy_dzeta*( dphi_dxi[i][p].slice(1) -
239  dphi_deta[i][p].slice(0) ) );
240 
241  curl_phi[i][p].slice(2) = inv_jac*( dz_dxi*( dphi_deta[i][p].slice(2) -
242  dphi_dzeta[i][p].slice(1) ) +
243  dz_deta*( dphi_dzeta[i][p].slice(0) -
244  dphi_dxi[i][p].slice(2) ) +
245  dz_dzeta*( dphi_dxi[i][p].slice(1) -
246  dphi_deta[i][p].slice(0) ) );
247  }
248 
249  break;
250  }
251 
252  default:
253  libmesh_error_msg("Invalid dim = " << dim);
254  } // switch(dim)
255 }
256 
257 template class LIBMESH_EXPORT HCurlFETransformation<RealGradient>;
258 
259 template<>
261 {
262  libmesh_error_msg("HCurl transformations only make sense for vector-valued elements.");
263 }
264 
265 template<>
267 {
268  libmesh_error_msg("HCurl transformations only make sense for vector-valued elements.");
269 }
270 
271 template<>
273 {
274  libmesh_error_msg("HCurl transformations only make sense for vector-valued elements.");
275 }
276 
277 template<>
278 void HCurlFETransformation<Real>::map_phi(const unsigned int,
279  const Elem * const,
280  const std::vector<Point> &,
281  const FEGenericBase<Real> &,
282  std::vector<std::vector<Real>> &,
283  bool) const
284 {
285  libmesh_error_msg("HCurl transformations only make sense for vector-valued elements.");
286 }
287 
288 template<>
290  const Elem * const,
291  const std::vector<Point> &,
292  const FEGenericBase<Real> &,
293  std::vector<std::vector<Real>> &) const
294 {
295  libmesh_error_msg("HCurl transformations only make sense for vector-valued elements.");
296 }
297 
298 
299 } // namespace libMesh
virtual void init_map_phi(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 for HCurl conforming...
unsigned int dim
const std::vector< Real > & get_detadz() const
Definition: fe_map.h:349
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
const std::vector< RealGradient > & get_dxyzdzeta() const
Definition: fe_map.h:255
The libMesh namespace provides an interface to certain functionality in the library.
const std::vector< std::vector< OutputShape > > & get_dphideta() const
Definition: fe_base.h:301
const std::vector< Real > & get_dxidz() const
Definition: fe_map.h:325
const std::vector< std::vector< OutputShape > > & get_dphidzeta() const
Definition: fe_base.h:309
FEType get_fe_type() const
Definition: fe_abstract.h:499
const std::vector< RealGradient > & get_dxyzdxi() const
Definition: fe_map.h:239
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
const std::vector< Real > & get_dzetadx() const
Definition: fe_map.h:357
const std::vector< Real > & get_jacobian() const
Definition: fe_map.h:223
const std::vector< Real > & get_dxidx() const
Definition: fe_map.h:309
const std::vector< Real > & get_dzetady() const
Definition: fe_map.h:365
const std::vector< Real > & get_dxidy() const
Definition: fe_map.h:317
const std::vector< Real > & get_dzetadz() const
Definition: fe_map.h:373
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 override
Evaluates the shape function curl in physical coordinates based on conforming finite element transfo...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< Real > & get_detady() const
Definition: fe_map.h:341
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 init_map_d2phi(const FEGenericBase< OutputShape > &fe) const override
Pre-requests any necessary data from FEMap.
const std::vector< RealGradient > & get_dxyzdeta() const
Definition: fe_map.h:247
const std::vector< std::vector< OutputShape > > & get_dphidxi() const
Definition: fe_base.h:293
const FEMap & get_fe_map() const
Definition: fe_abstract.h:533
const std::vector< std::vector< Real > > & get_d2xidxyz2() const
Second derivatives of "xi" reference coordinate wrt physical coordinates.
Definition: fe_map.h:381
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:111
This class forms the foundation from which generic finite elements may be derived.
virtual void init_map_dphi(const FEGenericBase< OutputShape > &fe) const override
Pre-requests any necessary data from FEMap.
const std::vector< Real > & get_detadx() const
Definition: fe_map.h:333