libMesh
fe_xyz.C
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 
19 
20 // Local includes
21 #include "libmesh/elem.h"
22 #include "libmesh/enum_to_string.h"
23 #include "libmesh/fe.h"
24 #include "libmesh/fe_interface.h"
25 #include "libmesh/fe_macro.h"
26 #include "libmesh/int_range.h"
27 #include "libmesh/libmesh_logging.h"
28 
29 namespace libMesh
30 {
31 
32 // ------------------------------------------------------------
33 // XYZ-specific implementations
34 
35 // Anonymous namespace for local helper functions
36 namespace {
37 
38 void xyz_nodal_soln(const Elem * elem,
39  const Order order,
40  const std::vector<Number> & elem_soln,
41  std::vector<Number> & nodal_soln,
42  const bool add_p_level)
43 {
44  const unsigned int n_nodes = elem->n_nodes();
45 
46  nodal_soln.resize(n_nodes);
47 
48  const Order totalorder = order + add_p_level*elem->p_level();
49 
50  switch (totalorder)
51  {
52  // Constant shape functions
53  case CONSTANT:
54  {
55  libmesh_assert_equal_to (elem_soln.size(), 1);
56 
57  std::fill(nodal_soln.begin(), nodal_soln.end(), elem_soln[0]);
58 
59  return;
60  }
61 
62 
63  // For other orders do interpolation at the nodes
64  // explicitly.
65  default:
66  {
67  // FEType object to be passed to various FEInterface functions below.
68  FEType fe_type(order, XYZ);
69 
70  const unsigned int n_sf =
71  FEInterface::n_shape_functions(fe_type, elem);
72  libmesh_assert_equal_to (elem_soln.size(), n_sf);
73 
74  // Zero before summation
75  std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
76 
77  for (unsigned int n=0; n<n_nodes; n++)
78  // u_i = Sum (alpha_i phi_i)
79  for (unsigned int i=0; i<n_sf; i++)
80  nodal_soln[n] += elem_soln[i] *
81  FEInterface::shape(fe_type, elem, i, elem->point(n));
82 
83  return;
84  } // default
85  } // switch
86 } // xyz_nodal_soln()
87 
88 
89 } // anonymous namespace
90 
91 
92 
93 
94 
95 
96 
97 template <unsigned int Dim>
98 void FEXYZ<Dim>::init_shape_functions(const std::vector<Point> & qp,
99  const Elem * elem)
100 {
101  libmesh_assert(elem);
102 
103  // FIXME: Is this redundant here? Who's calling init_shape_functions
104  // from code that hasn't already done a determine_calculations()?
105  this->determine_calculations();
106 
107  // Start logging the shape function initialization
108  LOG_SCOPE("init_shape_functions()", "FE");
109 
110  // The number of quadrature points.
111  const std::size_t n_qp = qp.size();
112 
113  // Number of shape functions in the finite element approximation
114  // space.
115  const unsigned int n_approx_shape_functions =
116  this->n_dofs(elem,
117  this->get_order());
118 
119  // resize the vectors to hold current data
120  // Phi are the shape functions used for the FE approximation
121  {
122  // (note: GCC 3.4.0 requires the use of this-> here)
123  if (this->calculate_phi)
124  this->phi.resize (n_approx_shape_functions);
125  if (this->calculate_dphi)
126  {
127  this->dphi.resize (n_approx_shape_functions);
128  this->dphidx.resize (n_approx_shape_functions);
129  this->dphidy.resize (n_approx_shape_functions);
130  this->dphidz.resize (n_approx_shape_functions);
131  }
132 
133 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
134  if (this->calculate_d2phi)
135  {
136  this->d2phi.resize (n_approx_shape_functions);
137  this->d2phidx2.resize (n_approx_shape_functions);
138  this->d2phidxdy.resize (n_approx_shape_functions);
139  this->d2phidxdz.resize (n_approx_shape_functions);
140  this->d2phidy2.resize (n_approx_shape_functions);
141  this->d2phidydz.resize (n_approx_shape_functions);
142  this->d2phidz2.resize (n_approx_shape_functions);
143  this->d2phidxi2.resize (n_approx_shape_functions);
144  }
145 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
146 
147  for (unsigned int i=0; i<n_approx_shape_functions; i++)
148  {
149  if (this->calculate_phi)
150  this->phi[i].resize (n_qp);
151  if (this->calculate_dphi)
152  {
153  this->dphi[i].resize (n_qp);
154  this->dphidx[i].resize (n_qp);
155  this->dphidy[i].resize (n_qp);
156  this->dphidz[i].resize (n_qp);
157  }
158 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
159  if (this->calculate_d2phi)
160  {
161  this->d2phi[i].resize (n_qp);
162  this->d2phidx2[i].resize (n_qp);
163  this->d2phidxdy[i].resize (n_qp);
164  this->d2phidxdz[i].resize (n_qp);
165  this->d2phidy2[i].resize (n_qp);
166  this->d2phidydz[i].resize (n_qp);
167  this->d2phidz2[i].resize (n_qp);
168  }
169 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
170  }
171  }
172 
173 
174 
175 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
176  //------------------------------------------------------------
177  // Initialize the data fields, which should only be used for infinite
178  // elements, to some sensible values, so that using a FE with the
179  // variational formulation of an InfFE, correct element matrices are
180  // returned
181 
182  {
183  this->weight.resize (n_qp);
184  this->dweight.resize (n_qp);
185  this->dphase.resize (n_qp);
186 
187  for (unsigned int p=0; p<n_qp; p++)
188  {
189  this->weight[p] = 1.;
190  this->dweight[p].zero();
191  this->dphase[p].zero();
192  }
193 
194  }
195 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
196 
197  if (this->calculate_dual)
198  this->init_dual_shape_functions(n_approx_shape_functions, n_qp);
199 }
200 
201 
202 
203 
204 template <unsigned int Dim>
206  const std::vector<Point> &)
207 {
208  libmesh_assert(elem);
209 
210  //-------------------------------------------------------------------------
211  // Compute the shape function values (and derivatives)
212  // at the Quadrature points. Note that the actual values
213  // have already been computed via init_shape_functions
214 
215  // Start logging the shape function computation
216  LOG_SCOPE("compute_shape_functions()", "FE");
217 
218  const std::vector<Point> & xyz_qp = this->get_xyz();
219 
220  // Compute the value of the derivative shape function i at quadrature point p
221  switch (this->dim)
222  {
223 
224  case 1:
225  {
226  if (this->calculate_phi)
227  for (auto i : index_range(this->phi))
228  for (auto p : index_range(this->phi[i]))
229  this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]);
230 
231  if (this->calculate_dphi)
232  for (auto i : index_range(this->dphi))
233  for (auto p : index_range(this->dphi[i]))
234  {
235  this->dphi[i][p](0) =
236  this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
237 
238  this->dphi[i][p](1) = this->dphidy[i][p] = 0.;
239  this->dphi[i][p](2) = this->dphidz[i][p] = 0.;
240  }
241 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
242  if (this->calculate_d2phi)
243  for (auto i : index_range(this->d2phi))
244  for (auto p : index_range(this->d2phi[i]))
245  {
246  this->d2phi[i][p](0,0) =
247  this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
248 
249 #if LIBMESH_DIM>1
250  this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] =
251  this->d2phi[i][p](1,0) = 0.;
252  this->d2phi[i][p](1,1) = this->d2phidy2[i][p] = 0.;
253 #if LIBMESH_DIM>2
254  this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] =
255  this->d2phi[i][p](2,0) = 0.;
256  this->d2phi[i][p](1,2) = this->d2phidydz[i][p] =
257  this->d2phi[i][p](2,1) = 0.;
258  this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = 0.;
259 #endif
260 #endif
261  }
262 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
263 
264  // All done
265  break;
266  }
267 
268  case 2:
269  {
270  if (this->calculate_phi)
271  for (auto i : index_range(this->phi))
272  for (auto p : index_range(this->phi[i]))
273  this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]);
274 
275  if (this->calculate_dphi)
276  for (auto i : index_range(this->dphi))
277  for (auto p : index_range(this->dphi[i]))
278  {
279  this->dphi[i][p](0) =
280  this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
281 
282  this->dphi[i][p](1) =
283  this->dphidy[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]);
284 
285 #if LIBMESH_DIM == 3
286  this->dphi[i][p](2) = // can only assign to the Z component if LIBMESH_DIM==3
287 #endif
288  this->dphidz[i][p] = 0.;
289  }
290 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
291  if (this->calculate_d2phi)
292  for (auto i : index_range(this->d2phi))
293  for (auto p : index_range(this->d2phi[i]))
294  {
295  this->d2phi[i][p](0,0) =
296  this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
297 
298  this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] =
299  this->d2phi[i][p](1,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]);
300  this->d2phi[i][p](1,1) =
301  this->d2phidy2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]);
302 #if LIBMESH_DIM>2
303  this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] =
304  this->d2phi[i][p](2,0) = 0.;
305  this->d2phi[i][p](1,2) = this->d2phidydz[i][p] =
306  this->d2phi[i][p](2,1) = 0.;
307  this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = 0.;
308 #endif
309  }
310 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
311 
312  // All done
313  break;
314  }
315 
316  case 3:
317  {
318  if (this->calculate_phi)
319  for (auto i : index_range(this->phi))
320  for (auto p : index_range(this->phi[i]))
321  this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]);
322 
323  if (this->calculate_dphi)
324  for (auto i : index_range(this->dphi))
325  for (auto p : index_range(this->dphi[i]))
326  {
327  this->dphi[i][p](0) =
328  this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
329 
330  this->dphi[i][p](1) =
331  this->dphidy[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]);
332 
333  this->dphi[i][p](2) =
334  this->dphidz[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]);
335  }
336 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
337  if (this->calculate_d2phi)
338  for (auto i : index_range(this->d2phi))
339  for (auto p : index_range(this->d2phi[i]))
340  {
341  this->d2phi[i][p](0,0) =
342  this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
343 
344  this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] =
345  this->d2phi[i][p](1,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]);
346  this->d2phi[i][p](1,1) =
347  this->d2phidy2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]);
348  this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] =
349  this->d2phi[i][p](2,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 3, xyz_qp[p]);
350  this->d2phi[i][p](1,2) = this->d2phidydz[i][p] =
351  this->d2phi[i][p](2,1) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 4, xyz_qp[p]);
352  this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 5, xyz_qp[p]);
353  }
354 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
355 
356  // All done
357  break;
358  }
359 
360  default:
361  libmesh_error_msg("ERROR: Invalid dimension " << this->dim);
362  }
363 }
364 
365 
366 // Instantiate (side_) nodal_soln() function for every dimension
367 LIBMESH_FE_NODAL_SOLN(XYZ, xyz_nodal_soln)
369 
370 
371 // Full specialization of n_dofs() function for every dimension
372 template <> unsigned int FE<0,XYZ>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
373 template <> unsigned int FE<1,XYZ>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
374 template <> unsigned int FE<2,XYZ>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
375 template <> unsigned int FE<3,XYZ>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
376 
377 template <> unsigned int FE<0,XYZ>::n_dofs(const Elem * e, const Order o) { return monomial_n_dofs(e, o); }
378 template <> unsigned int FE<1,XYZ>::n_dofs(const Elem * e, const Order o) { return monomial_n_dofs(e, o); }
379 template <> unsigned int FE<2,XYZ>::n_dofs(const Elem * e, const Order o) { return monomial_n_dofs(e, o); }
380 template <> unsigned int FE<3,XYZ>::n_dofs(const Elem * e, const Order o) { return monomial_n_dofs(e, o); }
381 
382 // Full specialization of n_dofs_at_node() function for every dimension.
383 // XYZ FEMs have no dofs at nodes, only element dofs.
384 template <> unsigned int FE<0,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
385 template <> unsigned int FE<1,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
386 template <> unsigned int FE<2,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
387 template <> unsigned int FE<3,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
388 
389 template <> unsigned int FE<0,XYZ>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
390 template <> unsigned int FE<1,XYZ>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
391 template <> unsigned int FE<2,XYZ>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
392 template <> unsigned int FE<3,XYZ>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
393 
394 // Full specialization of n_dofs_per_elem() function for every dimension.
395 template <> unsigned int FE<0,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
396 template <> unsigned int FE<1,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
397 template <> unsigned int FE<2,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
398 template <> unsigned int FE<3,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
399 
400 template <> unsigned int FE<0,XYZ>::n_dofs_per_elem(const Elem & e, const Order o) { return monomial_n_dofs(&e, o); }
401 template <> unsigned int FE<1,XYZ>::n_dofs_per_elem(const Elem & e, const Order o) { return monomial_n_dofs(&e, o); }
402 template <> unsigned int FE<2,XYZ>::n_dofs_per_elem(const Elem & e, const Order o) { return monomial_n_dofs(&e, o); }
403 template <> unsigned int FE<3,XYZ>::n_dofs_per_elem(const Elem & e, const Order o) { return monomial_n_dofs(&e, o); }
404 
405 // Full specialization of get_continuity() function for every dimension.
406 template <> FEContinuity FE<0,XYZ>::get_continuity() const { return DISCONTINUOUS; }
407 template <> FEContinuity FE<1,XYZ>::get_continuity() const { return DISCONTINUOUS; }
408 template <> FEContinuity FE<2,XYZ>::get_continuity() const { return DISCONTINUOUS; }
409 template <> FEContinuity FE<3,XYZ>::get_continuity() const { return DISCONTINUOUS; }
410 
411 // Full specialization of is_hierarchic() function for every dimension.
412 // The XYZ shape functions are hierarchic!
413 template <> bool FE<0,XYZ>::is_hierarchic() const { return true; }
414 template <> bool FE<1,XYZ>::is_hierarchic() const { return true; }
415 template <> bool FE<2,XYZ>::is_hierarchic() const { return true; }
416 template <> bool FE<3,XYZ>::is_hierarchic() const { return true; }
417 
418 #ifdef LIBMESH_ENABLE_AMR
419 
420 // Full specialization of compute_constraints() function for 2D and
421 // 3D only. There are no constraints for discontinuous elements, so
422 // we do nothing.
423 template <> void FE<2,XYZ>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem *) {}
424 template <> void FE<3,XYZ>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem *) {}
425 
426 #endif // #ifdef LIBMESH_ENABLE_AMR
427 
428 // Full specialization of shapes_need_reinit() function for every dimension.
429 template <> bool FE<0,XYZ>::shapes_need_reinit() const { return true; }
430 template <> bool FE<1,XYZ>::shapes_need_reinit() const { return true; }
431 template <> bool FE<2,XYZ>::shapes_need_reinit() const { return true; }
432 template <> bool FE<3,XYZ>::shapes_need_reinit() const { return true; }
433 
434 
435 // Explicit instantiations for non-static FEXYZ member functions.
436 // These non-static member functions map more naturally to explicit
437 // instantiations than the functions above:
438 //
439 // 1.) Since they are member functions, they rely on
440 // private/protected member data, and therefore do not work well
441 // with the "anonymous function call" model we've used above for
442 // the specializations.
443 //
444 // 2.) There is (IMHO) less chance of the linker calling the
445 // wrong version of one of these member functions, since there is
446 // only one FEXYZ.
447 template LIBMESH_EXPORT void FEXYZ<0>::init_shape_functions(const std::vector<Point> &, const Elem *);
448 template LIBMESH_EXPORT void FEXYZ<1>::init_shape_functions(const std::vector<Point> &, const Elem *);
449 template LIBMESH_EXPORT void FEXYZ<2>::init_shape_functions(const std::vector<Point> &, const Elem *);
450 template LIBMESH_EXPORT void FEXYZ<3>::init_shape_functions(const std::vector<Point> &, const Elem *);
451 
452 template LIBMESH_EXPORT void FEXYZ<0>::compute_shape_functions(const Elem *,const std::vector<Point> &);
453 template LIBMESH_EXPORT void FEXYZ<1>::compute_shape_functions(const Elem *,const std::vector<Point> &);
454 template LIBMESH_EXPORT void FEXYZ<2>::compute_shape_functions(const Elem *,const std::vector<Point> &);
455 template LIBMESH_EXPORT void FEXYZ<3>::compute_shape_functions(const Elem *,const std::vector<Point> &);
456 
457 } // namespace libMesh
static unsigned int n_dofs(const ElemType t, const Order o)
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
unsigned int dim
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
virtual bool shapes_need_reinit() const override
The libMesh namespace provides an interface to certain functionality in the library.
LIBMESH_FE_NODAL_SOLN(LIBMESH_FE_SIDE_NODAL_SOLN() LIBMESH_DEFAULT_NDOFS(BERNSTEIN) template<> FEContinuity FE< 0 BERNSTEIN, bernstein_nodal_soln)
Definition: fe_bernstein.C:415
virtual void compute_shape_functions(const Elem *elem, const std::vector< Point > &qp) override
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
Definition: fe_xyz.C:205
virtual bool is_hierarchic() const override
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
A specific instantiation of the FEBase class.
Definition: fe.h:127
LIBMESH_FE_SIDE_NODAL_SOLN(HIERARCHIC_VEC)
const dof_id_type n_nodes
Definition: tecplot_io.C:67
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:272
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:437
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
libmesh_assert(ctx)
unsigned int monomial_n_dofs(const ElemType t, const Order o)
Helper functions for Discontinuous-Pn type basis functions.
Definition: fe_monomial.C:37
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
virtual FEContinuity get_continuity() const override
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
virtual void init_shape_functions(const std::vector< Point > &qp, const Elem *e) override
Update the various member data fields phi, dphidxi, dphideta, dphidzeta, etc.
Definition: fe_xyz.C:98
The constraint matrix storage format.
Definition: dof_map.h:108
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
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)