libMesh
fe_boundary.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 
19 
20 // C++ includes
21 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
22 #include <cmath> // for std::sqrt
23 
24 
25 // Local includes
26 #include "libmesh/libmesh_common.h"
27 #include "libmesh/fe.h"
28 #include "libmesh/fe_interface.h"
29 #include "libmesh/quadrature.h"
30 #include "libmesh/elem.h"
31 #include "libmesh/libmesh_logging.h"
32 #include "libmesh/tensor_value.h" // May be necessary if destructors
33 // get instantiated here
34 
35 namespace libMesh
36 {
37 
38 //-------------------------------------------------------
39 // Full specializations for useless methods in 0D, 1D
40 #define REINIT_ERROR(_dim, _type, _func) \
41  template <> \
42  void FE<_dim,_type>::_func(const Elem *, \
43  const unsigned int, \
44  const Real, \
45  const std::vector<Point> * const, \
46  const std::vector<Real> * const)
47 
48 #define SIDEMAP_ERROR(_dim, _type, _func) \
49  template <> \
50  void FE<_dim,_type>::_func(const Elem *, \
51  const Elem *, \
52  const unsigned int, \
53  const std::vector<Point> &, \
54  std::vector<Point> &)
55 
56 #define FACE_EDGE_SHAPE_ERROR(_dim, _func) \
57  template <> \
58  void FEMap::_func<_dim>(const std::vector<Point> &, \
59  const Elem *) \
60  { \
61  libmesh_error_msg("ERROR: This method makes no sense for low-D elements!"); \
62  }
63 
64 // 0D error instantiations
65 #define LIBMESH_ERRORS_IN_LOW_D(_type) \
66 REINIT_ERROR(0, _type, reinit) { libmesh_error_msg("ERROR: Cannot reinit 0D " #_type " elements!"); } \
67 REINIT_ERROR(0, _type, edge_reinit) { libmesh_error_msg("ERROR: Cannot edge_reinit 0D " #_type " elements!"); } \
68 SIDEMAP_ERROR(0, _type, side_map) { libmesh_error_msg("ERROR: Cannot side_map 0D " #_type " elements!"); } \
69 SIDEMAP_ERROR(0, _type, edge_map) { libmesh_error_msg("ERROR: Cannot edge_map 0D " #_type " elements!"); } \
70 REINIT_ERROR(1, _type, edge_reinit) { libmesh_error_msg("ERROR: Cannot edge_reinit 1D " #_type " elements!"); } \
71 SIDEMAP_ERROR(1, _type, edge_map) { libmesh_error_msg("ERROR: Cannot edge_map 1D " #_type " elements!"); }
72 
91 
92 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
96 #endif
97 
98 // Special error instantiations
99 REINIT_ERROR(1, NEDELEC_ONE, reinit) { libmesh_error_msg("ERROR: Cannot reinit 1D NEDELEC_ONE elements!"); }
100 SIDEMAP_ERROR(1, NEDELEC_ONE, side_map) { libmesh_error_msg("ERROR: Cannot side_map 1D NEDELEC_ONE elements!"); }
101 REINIT_ERROR(1, RAVIART_THOMAS, reinit) { libmesh_error_msg("ERROR: Cannot reinit 1D RAVIART_THOMAS elements!"); }
102 SIDEMAP_ERROR(1, RAVIART_THOMAS, side_map) { libmesh_error_msg("ERROR: Cannot side_map 1D RAVIART_THOMAS elements!"); }
103 REINIT_ERROR(1, L2_RAVIART_THOMAS, reinit) { libmesh_error_msg("ERROR: Cannot reinit 1D L2_RAVIART_THOMAS elements!"); }
104 SIDEMAP_ERROR(1, L2_RAVIART_THOMAS, side_map) { libmesh_error_msg("ERROR: Cannot side_map 1D L2_RAVIART_THOMAS elements!"); }
105 
106 //-------------------------------------------------------
107 // Methods for 2D, 3D
108 template <unsigned int Dim, FEFamily T>
109 void FE<Dim,T>::reinit(const Elem * elem,
110  const unsigned int s,
111  const Real /* tolerance */,
112  const std::vector<Point> * const pts,
113  const std::vector<Real> * const weights)
114 {
115  libmesh_assert(elem);
116  libmesh_assert (this->qrule != nullptr || pts != nullptr);
117  // We now do this for 1D elements!
118  // libmesh_assert_not_equal_to (Dim, 1);
119 
120  // If we called this function redundantly (e.g. in an FEMContext
121  // that is asked not to do any side calculations) then let's skip the
122  // whole inverse_map process that calculates side points
123  if (this->calculating_nothing())
124  {
125  this->calculations_started = true; // Ironic
126  return;
127  }
128 
129  // We're (possibly re-) calculating now! FIXME - we currently
130  // expect to be able to use side_map and JxW later, but we could
131  // optimize further here.
132  this->_fe_map->add_calculations();
133  this->_fe_map->get_JxW();
134  this->_fe_map->get_xyz();
135  this->determine_calculations();
136 
137  // Build the side of interest
138  const std::unique_ptr<const Elem> side(elem->build_side_ptr(s));
139 
140  // Find the max p_level to select
141  // the right quadrature rule for side integration
142  unsigned int side_p_level = elem->p_level();
143  if (elem->neighbor_ptr(s) != nullptr)
144  side_p_level = std::max(side_p_level, elem->neighbor_ptr(s)->p_level());
145 
146  // Initialize the shape functions at the user-specified
147  // points
148  if (pts != nullptr)
149  {
150  // The shape functions do not correspond to the qrule
151  this->shapes_on_quadrature = false;
152 
153  // Initialize the face shape functions
154  this->_fe_map->template init_face_shape_functions<Dim>(*pts, side.get());
155 
156  // Compute the Jacobian*Weight on the face for integration
157  if (weights != nullptr)
158  {
159  this->_fe_map->compute_face_map (Dim, *weights, side.get());
160  }
161  else
162  {
163  std::vector<Real> dummy_weights (pts->size(), 1.);
164  this->_fe_map->compute_face_map (Dim, dummy_weights, side.get());
165  }
166  }
167  // If there are no user specified points, we use the
168  // quadrature rule
169  else
170  {
171  // initialize quadrature rule
172  this->qrule->init(side->type(), side_p_level);
173 
174  if (this->qrule->shapes_need_reinit())
175  this->shapes_on_quadrature = false;
176 
177  // FIXME - could this break if the same FE object was used
178  // for both volume and face integrals? - RHS
179  // We might not need to reinitialize the shape functions
180  if ((this->get_type() != elem->type()) ||
181  (side->type() != last_side) ||
182  (this->_elem_p_level != side_p_level) ||
183  this->shapes_need_reinit() ||
184  !this->shapes_on_quadrature)
185  {
186  // Set the element type and p_level
187  this->elem_type = elem->type();
188  this->_elem_p_level = side_p_level;
189 
190  // Set the last_side
191  last_side = side->type();
192 
193  // Set the last p level
194  this->_p_level = this->_add_p_level_in_reinit * side_p_level;
195 
196  // Initialize the face shape functions
197  this->_fe_map->template init_face_shape_functions<Dim>(this->qrule->get_points(), side.get());
198  }
199 
200  // Compute the Jacobian*Weight on the face for integration
201  this->_fe_map->compute_face_map (Dim, this->qrule->get_weights(), side.get());
202 
203  // The shape functions correspond to the qrule
204  this->shapes_on_quadrature = true;
205  }
206 
207  // make a copy of the Jacobian for integration
208  const std::vector<Real> JxW_int(this->_fe_map->get_JxW());
209 
210  // make a copy of shape on quadrature info
211  bool shapes_on_quadrature_side = this->shapes_on_quadrature;
212 
213  // Find where the integration points are located on the
214  // full element.
215  const std::vector<Point> * ref_qp;
216  if (pts != nullptr)
217  ref_qp = pts;
218  else
219  ref_qp = &this->qrule->get_points();
220 
221  std::vector<Point> qp;
222  this->side_map(elem, side.get(), s, *ref_qp, qp);
223 
224  // compute the shape function and derivative values
225  // at the points qp
226  this->reinit (elem, &qp);
227 
228  this->shapes_on_quadrature = shapes_on_quadrature_side;
229 
230  // copy back old data
231  this->_fe_map->get_JxW() = JxW_int;
232 }
233 
234 
235 
236 template <unsigned int Dim, FEFamily T>
237 void FE<Dim,T>::edge_reinit(const Elem * elem,
238  const unsigned int e,
239  const Real /* tolerance */,
240  const std::vector<Point> * const pts,
241  const std::vector<Real> * const weights)
242 {
243  libmesh_assert(elem);
244  libmesh_assert (this->qrule != nullptr || pts != nullptr);
245  // We don't do this for 1D elements!
246  libmesh_assert_not_equal_to (Dim, 1);
247 
248  // We're (possibly re-) calculating now! Time to determine what.
249  // FIXME - we currently just assume that we're using JxW and calling
250  // edge_map later.
251  this->_fe_map->add_calculations();
252  this->_fe_map->get_JxW();
253  this->_fe_map->get_xyz();
254  this->determine_calculations();
255 
256  // Build the side of interest
257  const std::unique_ptr<const Elem> edge(elem->build_edge_ptr(e));
258 
259  // Initialize the shape functions at the user-specified
260  // points
261  if (pts != nullptr)
262  {
263  // The shape functions do not correspond to the qrule
264  this->shapes_on_quadrature = false;
265 
266  // Initialize the edge shape functions
267  this->_fe_map->template init_edge_shape_functions<Dim> (*pts, edge.get());
268 
269  // Compute the Jacobian*Weight on the face for integration
270  if (weights != nullptr)
271  {
272  this->_fe_map->compute_edge_map (Dim, *weights, edge.get());
273  }
274  else
275  {
276  std::vector<Real> dummy_weights (pts->size(), 1.);
277  this->_fe_map->compute_edge_map (Dim, dummy_weights, edge.get());
278  }
279  }
280  // If there are no user specified points, we use the
281  // quadrature rule
282  else
283  {
284  // initialize quadrature rule
285  this->qrule->init(edge->type(), elem->p_level());
286 
287  if (this->qrule->shapes_need_reinit())
288  this->shapes_on_quadrature = false;
289 
290  // We might not need to reinitialize the shape functions
291  if ((this->get_type() != elem->type()) ||
292  (edge->type() != last_edge) ||
293  this->shapes_need_reinit() ||
294  !this->shapes_on_quadrature)
295  {
296  // Set the element type
297  this->elem_type = elem->type();
298 
299  // Set the last_edge
300  last_edge = edge->type();
301 
302  // Initialize the edge shape functions
303  this->_fe_map->template init_edge_shape_functions<Dim> (this->qrule->get_points(), edge.get());
304  }
305 
306  // Compute the Jacobian*Weight on the face for integration
307  this->_fe_map->compute_edge_map (Dim, this->qrule->get_weights(), edge.get());
308 
309  // The shape functions correspond to the qrule
310  this->shapes_on_quadrature = true;
311  }
312 
313  // make a copy of the Jacobian for integration
314  const std::vector<Real> JxW_int(this->_fe_map->get_JxW());
315 
316  // Find where the integration points are located on the
317  // full element.
318  const std::vector<Point> * ref_qp;
319  if (pts != nullptr)
320  ref_qp = pts;
321  else
322  ref_qp = & this->qrule->get_points();
323 
324  std::vector<Point> qp;
325  this->edge_map(elem, edge.get(), e, *ref_qp, qp);
326 
327  // compute the shape function and derivative values
328  // at the points qp
329  this->reinit (elem, &qp);
330 
331  // copy back old data
332  this->_fe_map->get_JxW() = JxW_int;
333 }
334 
335 template <unsigned int Dim, FEFamily T>
336 void FE<Dim,T>::side_map (const Elem * elem,
337  const Elem * side,
338  const unsigned int s,
339  const std::vector<Point> & reference_side_points,
340  std::vector<Point> & reference_points)
341 {
342  // We're calculating mappings - we need at least first order info
343  this->calculate_phi = true;
344  this->determine_calculations();
345 
346  unsigned int side_p_level = elem->p_level();
347  if (elem->neighbor_ptr(s) != nullptr)
348  side_p_level = std::max(side_p_level, elem->neighbor_ptr(s)->p_level());
349 
350  if (side->type() != last_side ||
351  side_p_level != this->_elem_p_level ||
352  !this->shapes_on_quadrature)
353  {
354  // Set the element type
355  this->elem_type = elem->type();
356  this->_elem_p_level = side_p_level;
357  this->_p_level = this->_add_p_level_in_reinit * side_p_level;
358 
359  // Set the last_side
360  last_side = side->type();
361 
362  // Initialize the face shape functions
363  this->_fe_map->template init_face_shape_functions<Dim>(reference_side_points, side);
364  }
365 
366  const unsigned int n_points =
367  cast_int<unsigned int>(reference_side_points.size());
368  reference_points.resize(n_points);
369  for (unsigned int i = 0; i < n_points; i++)
370  reference_points[i].zero();
371 
372  std::vector<Point> refspace_nodes;
373  this->get_refspace_nodes(elem->type(), refspace_nodes);
374 
375  const std::vector<std::vector<Real>> & psi_map = this->_fe_map->get_psi();
376 
377  // sum over the nodes
378  for (auto i : index_range(psi_map))
379  {
380  const Point & side_node = refspace_nodes[elem->local_side_node(s,i)];
381  for (unsigned int p=0; p<n_points; p++)
382  reference_points[p].add_scaled (side_node, psi_map[i][p]);
383  }
384 }
385 
386 template <unsigned int Dim, FEFamily T>
387 void FE<Dim,T>::edge_map (const Elem * elem,
388  const Elem * edge,
389  const unsigned int e,
390  const std::vector<Point> & reference_edge_points,
391  std::vector<Point> & reference_points)
392 {
393  // We're calculating mappings - we need at least first order info
394  this->calculate_phi = true;
395  this->determine_calculations();
396 
397  unsigned int edge_p_level = elem->p_level();
398 
399  if (edge->type() != last_edge ||
400  edge_p_level != this->_elem_p_level ||
401  !this->shapes_on_quadrature)
402  {
403  // Set the element type
404  this->elem_type = elem->type();
405  this->_elem_p_level = edge_p_level;
406  this->_p_level = this->_add_p_level_in_reinit * edge_p_level;
407 
408  // Set the last_edge
409  last_edge = edge->type();
410 
411  // Initialize the edge shape functions
412  this->_fe_map->template init_edge_shape_functions<Dim>(reference_edge_points, edge);
413  }
414 
415  const unsigned int n_points =
416  cast_int<unsigned int>(reference_edge_points.size());
417  reference_points.resize(n_points);
418  for (unsigned int i = 0; i < n_points; i++)
419  reference_points[i].zero();
420 
421  std::vector<Point> refspace_nodes;
422  this->get_refspace_nodes(elem->type(), refspace_nodes);
423 
424  const std::vector<std::vector<Real>> & psi_map = this->_fe_map->get_psi();
425 
426  // sum over the nodes
427  for (auto i : index_range(psi_map))
428  {
429  const Point & edge_node = refspace_nodes[elem->local_edge_node(e,i)];
430  for (unsigned int p=0; p<n_points; p++)
431  reference_points[p].add_scaled (edge_node, psi_map[i][p]);
432  }
433 }
434 
435 
436 template<unsigned int Dim>
437 void FEMap::init_face_shape_functions(const std::vector<Point> & qp,
438  const Elem * side)
439 {
440  // Start logging the shape function initialization
441  LOG_SCOPE("init_face_shape_functions()", "FEMap");
442 
443  libmesh_assert(side);
444 
445  // We're calculating now!
446  this->determine_calculations();
447 
448  // The element type and order to use in
449  // the map
450  const FEFamily mapping_family = FEMap::map_fe_type(*side);
451  const FEType map_fe_type(side->default_order(), mapping_family);
452 
453  // The number of quadrature points.
454  const unsigned int n_qp = cast_int<unsigned int>(qp.size());
455 
456  // Do not use the p_level(), if any, that is inherited by the side.
457  const unsigned int n_mapping_shape_functions =
458  FEInterface::n_shape_functions(map_fe_type, /*extra_order=*/0, side);
459 
460  // resize the vectors to hold current data
461  // Psi are the shape functions used for the FE mapping
462  if (calculate_xyz)
463  this->psi_map.resize (n_mapping_shape_functions);
464 
465  if (Dim > 1)
466  {
467  if (calculate_dxyz)
468  this->dpsidxi_map.resize (n_mapping_shape_functions);
469 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
470  if (calculate_d2xyz)
471  this->d2psidxi2_map.resize (n_mapping_shape_functions);
472 #endif
473  }
474 
475  if (Dim == 3)
476  {
477  if (calculate_dxyz)
478  this->dpsideta_map.resize (n_mapping_shape_functions);
479 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
480  if (calculate_d2xyz)
481  {
482  this->d2psidxideta_map.resize (n_mapping_shape_functions);
483  this->d2psideta2_map.resize (n_mapping_shape_functions);
484  }
485 #endif
486  }
487 
488  FEInterface::shape_ptr shape_ptr =
490 
491  FEInterface::shape_deriv_ptr shape_deriv_ptr =
493 
494 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
495  FEInterface::shape_second_deriv_ptr shape_second_deriv_ptr =
497 #endif
498 
499  for (unsigned int i=0; i<n_mapping_shape_functions; i++)
500  {
501  // Allocate space to store the values of the shape functions
502  // and their first and second derivatives at the quadrature points.
503  if (calculate_xyz)
504  this->psi_map[i].resize (n_qp);
505  if (Dim > 1)
506  {
507  if (calculate_dxyz)
508  this->dpsidxi_map[i].resize (n_qp);
509 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
510  if (calculate_d2xyz)
511  this->d2psidxi2_map[i].resize (n_qp);
512 #endif
513  }
514  if (Dim == 3)
515  {
516  if (calculate_dxyz)
517  this->dpsideta_map[i].resize (n_qp);
518 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
519  if (calculate_d2xyz)
520  {
521  this->d2psidxideta_map[i].resize (n_qp);
522  this->d2psideta2_map[i].resize (n_qp);
523  }
524 #endif
525  }
526 
527 
528  // Compute the value of mapping shape function i, and its first
529  // and second derivatives at quadrature point p
530  for (unsigned int p=0; p<n_qp; p++)
531  {
532  if (calculate_xyz)
533  this->psi_map[i][p] = shape_ptr (map_fe_type, side, i, qp[p], false);
534  if (Dim > 1)
535  {
536  if (calculate_dxyz)
537  this->dpsidxi_map[i][p] = shape_deriv_ptr (map_fe_type, side, i, 0, qp[p], false);
538 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
539  if (calculate_d2xyz)
540  this->d2psidxi2_map[i][p] = shape_second_deriv_ptr(map_fe_type, side, i, 0, qp[p], false);
541 #endif
542  }
543  // libMesh::out << "this->d2psidxi2_map["<<i<<"][p]=" << d2psidxi2_map[i][p] << std::endl;
544 
545  // If we are in 3D, then our sides are 2D faces.
546  // For the second derivatives, we must also compute the cross
547  // derivative d^2() / dxi deta
548  if (Dim == 3)
549  {
550  if (calculate_dxyz)
551  this->dpsideta_map[i][p] = shape_deriv_ptr (map_fe_type, side, i, 1, qp[p], false);
552 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
553  if (calculate_d2xyz)
554  {
555  this->d2psidxideta_map[i][p] = shape_second_deriv_ptr(map_fe_type, side, i, 1, qp[p], false);
556  this->d2psideta2_map[i][p] = shape_second_deriv_ptr(map_fe_type, side, i, 2, qp[p], false);
557  }
558 #endif
559  }
560  }
561  }
562 }
563 
564 template<unsigned int Dim>
565 void FEMap::init_edge_shape_functions(const std::vector<Point> & qp,
566  const Elem * edge)
567 {
568  // Start logging the shape function initialization
569  LOG_SCOPE("init_edge_shape_functions()", "FEMap");
570 
571  libmesh_assert(edge);
572 
573  // We're calculating now!
574  this->determine_calculations();
575 
576  // The element type and order to use in
577  // the map
578  const FEFamily mapping_family = FEMap::map_fe_type(*edge);
579  const FEType map_fe_type(edge->default_order(), mapping_family);
580 
581  // The number of quadrature points.
582  const unsigned int n_qp = cast_int<unsigned int>(qp.size());
583 
584  // Do not use the p_level(), if any, that is inherited by the side.
585  const unsigned int n_mapping_shape_functions =
586  FEInterface::n_shape_functions(map_fe_type, /*extra_order=*/0, edge);
587 
588  // resize the vectors to hold current data
589  // Psi are the shape functions used for the FE mapping
590  if (calculate_xyz)
591  this->psi_map.resize (n_mapping_shape_functions);
592  if (calculate_dxyz)
593  this->dpsidxi_map.resize (n_mapping_shape_functions);
594 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
595  if (calculate_d2xyz)
596  this->d2psidxi2_map.resize (n_mapping_shape_functions);
597 #endif
598 
599  FEInterface::shape_ptr shape_ptr =
601 
602  FEInterface::shape_deriv_ptr shape_deriv_ptr =
604 
605 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
606  FEInterface::shape_second_deriv_ptr shape_second_deriv_ptr =
608 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
609 
610  for (unsigned int i=0; i<n_mapping_shape_functions; i++)
611  {
612  // Allocate space to store the values of the shape functions
613  // and their first and second derivatives at the quadrature points.
614  if (calculate_xyz)
615  this->psi_map[i].resize (n_qp);
616  if (calculate_dxyz)
617  this->dpsidxi_map[i].resize (n_qp);
618 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
619  if (calculate_d2xyz)
620  this->d2psidxi2_map[i].resize (n_qp);
621 #endif
622 
623  // Compute the value of mapping shape function i, and its first
624  // and second derivatives at quadrature point p
625  for (unsigned int p=0; p<n_qp; p++)
626  {
627  if (calculate_xyz)
628  this->psi_map[i][p] = shape_ptr (map_fe_type, edge, i, qp[p], false);
629  if (calculate_dxyz)
630  this->dpsidxi_map[i][p] = shape_deriv_ptr (map_fe_type, edge, i, 0, qp[p], false);
631 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
632  if (calculate_d2xyz)
633  this->d2psidxi2_map[i][p] = shape_second_deriv_ptr(map_fe_type, edge, i, 0, qp[p], false);
634 #endif
635  }
636  }
637 }
638 
639 
640 
641 void FEMap::compute_face_map(int dim, const std::vector<Real> & qw,
642  const Elem * side)
643 {
644  libmesh_assert(side);
645 
646  // We're calculating now!
647  this->determine_calculations();
648 
649  LOG_SCOPE("compute_face_map()", "FEMap");
650 
651  // The number of quadrature points.
652  const unsigned int n_qp = cast_int<unsigned int>(qw.size());
653 
654  const FEFamily mapping_family = FEMap::map_fe_type(*side);
655  const Order mapping_order (side->default_order());
656  const FEType map_fe_type(mapping_order, mapping_family);
657 
658  // Do not use the p_level(), if any, that is inherited by the side.
659  const unsigned int n_mapping_shape_functions =
660  FEInterface::n_shape_functions(map_fe_type, /*extra_order=*/0, side);
661 
662  switch (dim)
663  {
664  case 1:
665  {
666  // A 1D finite element, currently assumed to be in 1D space
667  // This means the boundary is a "0D finite element", a
668  // NODEELEM.
669 
670  // Resize the vectors to hold data at the quadrature points
671  {
672  if (calculate_xyz)
673  this->xyz.resize(n_qp);
674  if (calculate_dxyz)
675  normals.resize(n_qp);
676 
677  if (calculate_dxyz)
678  this->JxW.resize(n_qp);
679  }
680 
681  // If we have no quadrature points, there's nothing else to do
682  if (!n_qp)
683  break;
684 
685  // We need to look back at the full edge to figure out the normal
686  // vector
687  const Elem * elem = side->interior_parent();
688  libmesh_assert (elem);
689  if (calculate_dxyz)
690  {
691  if (side->node_id(0) == elem->node_id(0))
692  normals[0] = Point(-1.);
693  else
694  {
695  libmesh_assert_equal_to (side->node_id(0),
696  elem->node_id(1));
697  normals[0] = Point(1.);
698  }
699  }
700 
701  // Calculate x at the point
702  if (calculate_xyz)
703  libmesh_assert_equal_to (this->psi_map.size(), 1);
704  // In the unlikely event we have multiple quadrature
705  // points, they'll be in the same place
706  for (unsigned int p=0; p<n_qp; p++)
707  {
708  if (calculate_xyz)
709  {
710  this->xyz[p].zero();
711  this->xyz[p].add_scaled (side->point(0), this->psi_map[0][p]);
712  }
713  if (calculate_dxyz)
714  {
715  normals[p] = normals[0];
716  this->JxW[p] = 1.0*qw[p];
717  }
718  }
719 
720  // done computing the map
721  break;
722  }
723 
724  case 2:
725  {
726  // A 2D finite element living in either 2D or 3D space.
727  // This means the boundary is a 1D finite element, i.e.
728  // and EDGE2 or EDGE3.
729  // Resize the vectors to hold data at the quadrature points
730  {
731  if (calculate_xyz)
732  this->xyz.resize(n_qp);
733  if (calculate_dxyz)
734  {
735  this->dxyzdxi_map.resize(n_qp);
736  this->tangents.resize(n_qp);
737  this->normals.resize(n_qp);
738 
739  this->JxW.resize(n_qp);
740  }
741 
742 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
743  if (calculate_d2xyz)
744  {
745  this->d2xyzdxi2_map.resize(n_qp);
746  this->curvatures.resize(n_qp);
747  }
748 #endif
749  }
750 
751  // Clear the entities that will be summed
752  // Compute the tangent & normal at the quadrature point
753  for (unsigned int p=0; p<n_qp; p++)
754  {
755  if (calculate_xyz)
756  this->xyz[p].zero();
757  if (calculate_dxyz)
758  {
759  this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D
760  this->dxyzdxi_map[p].zero();
761  }
762 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
763  if (calculate_d2xyz)
764  this->d2xyzdxi2_map[p].zero();
765 #endif
766  }
767 
768  // compute x, dxdxi at the quadrature points
769  for (unsigned int i=0; i<n_mapping_shape_functions; i++) // sum over the nodes
770  {
771  const Point & side_point = side->point(i);
772 
773  for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
774  {
775  if (calculate_xyz)
776  this->xyz[p].add_scaled (side_point, this->psi_map[i][p]);
777  if (calculate_dxyz)
778  this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]);
779 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
780  if (calculate_d2xyz)
781  this->d2xyzdxi2_map[p].add_scaled(side_point, this->d2psidxi2_map[i][p]);
782 #endif
783  }
784  }
785 
786  // Compute the tangent & normal at the quadrature point
787  if (calculate_dxyz)
788  {
789  for (unsigned int p=0; p<n_qp; p++)
790  {
791  // The first tangent comes from just the edge's Jacobian
792  this->tangents[p][0] = this->dxyzdxi_map[p].unit();
793 
794 #if LIBMESH_DIM == 2
795  // For a 2D element living in 2D, the normal is given directly
796  // from the entries in the edge Jacobian.
797  this->normals[p] = (Point(this->dxyzdxi_map[p](1), -this->dxyzdxi_map[p](0), 0.)).unit();
798 
799 #elif LIBMESH_DIM == 3
800  // For a 2D element living in 3D, there is a second tangent.
801  // For the second tangent, we need to refer to the full
802  // element's (not just the edge's) Jacobian.
803  const Elem * elem = side->interior_parent();
804  libmesh_assert(elem);
805 
806  // Inverse map xyz[p] to a reference point on the parent...
807  Point reference_point = FEMap::inverse_map(2, elem, this->xyz[p]);
808 
809  // Get dxyz/dxi and dxyz/deta from the parent map.
810  Point dx_dxi = FEMap::map_deriv (2, elem, 0, reference_point);
811  Point dx_deta = FEMap::map_deriv (2, elem, 1, reference_point);
812 
813  // The second tangent vector is formed by crossing these vectors.
814  tangents[p][1] = dx_dxi.cross(dx_deta).unit();
815 
816  // Finally, the normal in this case is given by crossing these
817  // two tangents.
818  normals[p] = tangents[p][0].cross(tangents[p][1]).unit();
819 #endif
820 
821 
822 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
823  // The curvature is computed via the familiar Frenet formula:
824  // curvature = [d^2(x) / d (xi)^2] dot [normal]
825  // For a reference, see:
826  // F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill, p. 310
827  //
828  // Note: The sign convention here is different from the
829  // 3D case. Concave-upward curves (smiles) have a positive
830  // curvature. Concave-downward curves (frowns) have a
831  // negative curvature. Be sure to take that into account!
832  if (calculate_d2xyz)
833  {
834  const Real numerator = this->d2xyzdxi2_map[p] * this->normals[p];
835  const Real denominator = this->dxyzdxi_map[p].norm_sq();
836  libmesh_assert_not_equal_to (denominator, 0);
837  curvatures[p] = numerator / denominator;
838  }
839 #endif
840  }
841 
842  // compute the jacobian at the quadrature points
843  for (unsigned int p=0; p<n_qp; p++)
844  {
845  const Real the_jac = this->dxyzdxi_map[p].norm();
846 
847  libmesh_assert_greater (the_jac, 0.);
848 
849  this->JxW[p] = the_jac*qw[p];
850  }
851  }
852 
853  // done computing the map
854  break;
855  }
856 
857 
858 
859  case 3:
860  {
861  // A 3D finite element living in 3D space.
862  // Resize the vectors to hold data at the quadrature points
863  {
864  if (calculate_xyz)
865  this->xyz.resize(n_qp);
866  if (calculate_dxyz)
867  {
868  this->dxyzdxi_map.resize(n_qp);
869  this->dxyzdeta_map.resize(n_qp);
870  this->tangents.resize(n_qp);
871  this->normals.resize(n_qp);
872  this->JxW.resize(n_qp);
873  }
874 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
875  if (calculate_d2xyz)
876  {
877  this->d2xyzdxi2_map.resize(n_qp);
878  this->d2xyzdxideta_map.resize(n_qp);
879  this->d2xyzdeta2_map.resize(n_qp);
880  this->curvatures.resize(n_qp);
881  }
882 #endif
883  }
884 
885  // Clear the entities that will be summed
886  for (unsigned int p=0; p<n_qp; p++)
887  {
888  if (calculate_xyz)
889  this->xyz[p].zero();
890  if (calculate_dxyz)
891  {
892  this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D
893  this->dxyzdxi_map[p].zero();
894  this->dxyzdeta_map[p].zero();
895  }
896 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
897  if (calculate_d2xyz)
898  {
899  this->d2xyzdxi2_map[p].zero();
900  this->d2xyzdxideta_map[p].zero();
901  this->d2xyzdeta2_map[p].zero();
902  }
903 #endif
904  }
905 
906  // compute x, dxdxi at the quadrature points
907  for (unsigned int i=0; i<n_mapping_shape_functions; i++) // sum over the nodes
908  {
909  const Point & side_point = side->point(i);
910 
911  for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
912  {
913  if (calculate_xyz)
914  this->xyz[p].add_scaled (side_point, this->psi_map[i][p]);
915  if (calculate_dxyz)
916  {
917  this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]);
918  this->dxyzdeta_map[p].add_scaled(side_point, this->dpsideta_map[i][p]);
919  }
920 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
921  if (calculate_d2xyz)
922  {
923  this->d2xyzdxi2_map[p].add_scaled (side_point, this->d2psidxi2_map[i][p]);
924  this->d2xyzdxideta_map[p].add_scaled(side_point, this->d2psidxideta_map[i][p]);
925  this->d2xyzdeta2_map[p].add_scaled (side_point, this->d2psideta2_map[i][p]);
926  }
927 #endif
928  }
929  }
930 
931  // Compute the tangents, normal, and curvature at the quadrature point
932  if (calculate_dxyz)
933  {
934  for (unsigned int p=0; p<n_qp; p++)
935  {
936  const Point n = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]);
937  this->normals[p] = n.unit();
938  this->tangents[p][0] = this->dxyzdxi_map[p].unit();
939  this->tangents[p][1] = n.cross(this->dxyzdxi_map[p]).unit();
940 
941 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
942  if (calculate_d2xyz)
943  {
944  // Compute curvature using the typical nomenclature
945  // of the first and second fundamental forms.
946  // For reference, see:
947  // 1) http://mathworld.wolfram.com/MeanCurvature.html
948  // (note -- they are using inward normal)
949  // 2) F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill
950  const Real L = -this->d2xyzdxi2_map[p] * this->normals[p];
951  const Real M = -this->d2xyzdxideta_map[p] * this->normals[p];
952  const Real N = -this->d2xyzdeta2_map[p] * this->normals[p];
953  const Real E = this->dxyzdxi_map[p].norm_sq();
954  const Real F = this->dxyzdxi_map[p] * this->dxyzdeta_map[p];
955  const Real G = this->dxyzdeta_map[p].norm_sq();
956 
957  const Real numerator = E*N -2.*F*M + G*L;
958  const Real denominator = E*G - F*F;
959  libmesh_assert_not_equal_to (denominator, 0.);
960  curvatures[p] = 0.5*numerator/denominator;
961  }
962 #endif
963  }
964 
965  // compute the jacobian at the quadrature points, see
966  // http://sp81.msi.umn.edu:999/fluent/fidap/help/theory/thtoc.htm
967  for (unsigned int p=0; p<n_qp; p++)
968  {
969  const Real g11 = (dxdxi_map(p)*dxdxi_map(p) +
970  dydxi_map(p)*dydxi_map(p) +
971  dzdxi_map(p)*dzdxi_map(p));
972 
973  const Real g12 = (dxdxi_map(p)*dxdeta_map(p) +
974  dydxi_map(p)*dydeta_map(p) +
975  dzdxi_map(p)*dzdeta_map(p));
976 
977  const Real g21 = g12;
978 
979  const Real g22 = (dxdeta_map(p)*dxdeta_map(p) +
980  dydeta_map(p)*dydeta_map(p) +
981  dzdeta_map(p)*dzdeta_map(p));
982 
983 
984  const Real the_jac = std::sqrt(g11*g22 - g12*g21);
985 
986  libmesh_assert_greater (the_jac, 0.);
987 
988  this->JxW[p] = the_jac*qw[p];
989  }
990  }
991 
992  // done computing the map
993  break;
994  }
995 
996 
997  default:
998  libmesh_error_msg("Invalid dimension dim = " << dim);
999  }
1000 }
1001 
1002 
1003 
1004 
1006  const std::vector<Real> & qw,
1007  const Elem * edge)
1008 {
1009  libmesh_assert(edge);
1010 
1011  if (dim == 2)
1012  {
1013  // A 2D finite element living in either 2D or 3D space.
1014  // The edges here are the sides of the element, so the
1015  // (misnamed) compute_face_map function does what we want
1016  this->compute_face_map(dim, qw, edge);
1017  return;
1018  }
1019 
1020  libmesh_assert_equal_to (dim, 3); // 1D is unnecessary and currently unsupported
1021 
1022  LOG_SCOPE("compute_edge_map()", "FEMap");
1023 
1024  // We're calculating now!
1025  this->determine_calculations();
1026 
1027  // The number of quadrature points.
1028  const unsigned int n_qp = cast_int<unsigned int>(qw.size());
1029 
1030  // Resize the vectors to hold data at the quadrature points
1031  if (calculate_xyz)
1032  this->xyz.resize(n_qp);
1033  if (calculate_dxyz)
1034  {
1035  this->dxyzdxi_map.resize(n_qp);
1036  this->dxyzdeta_map.resize(n_qp);
1037  this->tangents.resize(n_qp);
1038  this->normals.resize(n_qp);
1039  this->JxW.resize(n_qp);
1040  }
1041 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1042  if (calculate_d2xyz)
1043  {
1044  this->d2xyzdxi2_map.resize(n_qp);
1045  this->d2xyzdxideta_map.resize(n_qp);
1046  this->d2xyzdeta2_map.resize(n_qp);
1047  this->curvatures.resize(n_qp);
1048  }
1049 #endif
1050 
1051  // Clear the entities that will be summed
1052  for (unsigned int p=0; p<n_qp; p++)
1053  {
1054  if (calculate_xyz)
1055  this->xyz[p].zero();
1056  if (calculate_dxyz)
1057  {
1058  this->tangents[p].resize(1);
1059  this->dxyzdxi_map[p].zero();
1060  this->dxyzdeta_map[p].zero();
1061  }
1062 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1063  if (calculate_d2xyz)
1064  {
1065  this->d2xyzdxi2_map[p].zero();
1066  this->d2xyzdxideta_map[p].zero();
1067  this->d2xyzdeta2_map[p].zero();
1068  }
1069 #endif
1070  }
1071 
1072  // compute x, dxdxi at the quadrature points
1073  for (unsigned int i=0,
1074  psi_map_size=cast_int<unsigned int>(psi_map.size());
1075  i != psi_map_size; i++) // sum over the nodes
1076  {
1077  const Point & edge_point = edge->point(i);
1078 
1079  for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
1080  {
1081  if (calculate_xyz)
1082  this->xyz[p].add_scaled (edge_point, this->psi_map[i][p]);
1083  if (calculate_dxyz)
1084  this->dxyzdxi_map[p].add_scaled (edge_point, this->dpsidxi_map[i][p]);
1085 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1086  if (calculate_d2xyz)
1087  this->d2xyzdxi2_map[p].add_scaled (edge_point, this->d2psidxi2_map[i][p]);
1088 #endif
1089  }
1090  }
1091 
1092  // Compute the tangents at the quadrature point
1093  // FIXME: normals (plural!) and curvatures are uncalculated
1094  if (calculate_dxyz)
1095  for (unsigned int p=0; p<n_qp; p++)
1096  {
1097  // const Point n = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]);
1098  this->tangents[p][0] = this->dxyzdxi_map[p].unit();
1099 
1100  // compute the jacobian at the quadrature points
1101  const Real the_jac = std::sqrt(this->dxdxi_map(p)*this->dxdxi_map(p) +
1102  this->dydxi_map(p)*this->dydxi_map(p) +
1103  this->dzdxi_map(p)*this->dzdxi_map(p));
1104 
1105  libmesh_assert_greater (the_jac, 0.);
1106 
1107  this->JxW[p] = the_jac*qw[p];
1108  }
1109 }
1110 
1111 
1112 // Explicit FEMap Instantiations
1113 FACE_EDGE_SHAPE_ERROR(0,init_face_shape_functions)
1114 template LIBMESH_EXPORT void FEMap::init_face_shape_functions<1>(const std::vector<Point> &, const Elem *);
1115 template LIBMESH_EXPORT void FEMap::init_face_shape_functions<2>(const std::vector<Point> &, const Elem *);
1116 template LIBMESH_EXPORT void FEMap::init_face_shape_functions<3>(const std::vector<Point> &, const Elem *);
1117 
1118 FACE_EDGE_SHAPE_ERROR(0,init_edge_shape_functions)
1119 template LIBMESH_EXPORT void FEMap::init_edge_shape_functions<1>(const std::vector<Point> &, const Elem *);
1120 template LIBMESH_EXPORT void FEMap::init_edge_shape_functions<2>(const std::vector<Point> &, const Elem *);
1121 template LIBMESH_EXPORT void FEMap::init_edge_shape_functions<3>(const std::vector<Point> &, const Elem *);
1122 
1123 //--------------------------------------------------------------
1124 // Explicit FE instantiations
1125 #define REINIT_AND_SIDE_MAPS(_type) \
1126 template LIBMESH_EXPORT void FE<1,_type>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const); \
1127 template LIBMESH_EXPORT void FE<1,_type>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &); \
1128 template LIBMESH_EXPORT void FE<2,_type>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const); \
1129 template LIBMESH_EXPORT void FE<2,_type>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &); \
1130 template LIBMESH_EXPORT void FE<2,_type>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &); \
1131 template LIBMESH_EXPORT void FE<2,_type>::edge_reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const); \
1132 template LIBMESH_EXPORT void FE<3,_type>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const); \
1133 template LIBMESH_EXPORT void FE<3,_type>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &); \
1134 template LIBMESH_EXPORT void FE<3,_type>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &); \
1135 template LIBMESH_EXPORT void FE<3,_type>::edge_reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const)
1136 
1137 REINIT_AND_SIDE_MAPS(LAGRANGE);
1138 REINIT_AND_SIDE_MAPS(LAGRANGE_VEC);
1139 REINIT_AND_SIDE_MAPS(L2_LAGRANGE);
1140 REINIT_AND_SIDE_MAPS(L2_LAGRANGE_VEC);
1141 REINIT_AND_SIDE_MAPS(HIERARCHIC);
1142 REINIT_AND_SIDE_MAPS(HIERARCHIC_VEC);
1143 REINIT_AND_SIDE_MAPS(L2_HIERARCHIC);
1144 REINIT_AND_SIDE_MAPS(L2_HIERARCHIC_VEC);
1145 REINIT_AND_SIDE_MAPS(SIDE_HIERARCHIC);
1146 REINIT_AND_SIDE_MAPS(CLOUGH);
1147 REINIT_AND_SIDE_MAPS(HERMITE);
1148 REINIT_AND_SIDE_MAPS(MONOMIAL);
1149 REINIT_AND_SIDE_MAPS(MONOMIAL_VEC);
1150 REINIT_AND_SIDE_MAPS(SCALAR);
1151 REINIT_AND_SIDE_MAPS(XYZ);
1152 
1153 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
1154 REINIT_AND_SIDE_MAPS(BERNSTEIN);
1155 REINIT_AND_SIDE_MAPS(SZABAB);
1156 REINIT_AND_SIDE_MAPS(RATIONAL_BERNSTEIN);
1157 #endif
1158 
1159 template LIBMESH_EXPORT void FE<2,SUBDIVISION>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1160 template LIBMESH_EXPORT void FE<2,SUBDIVISION>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1161 
1162 template LIBMESH_EXPORT void FE<2,NEDELEC_ONE>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1163 template LIBMESH_EXPORT void FE<2,NEDELEC_ONE>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1164 template LIBMESH_EXPORT void FE<2,NEDELEC_ONE>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1165 template LIBMESH_EXPORT void FE<2,NEDELEC_ONE>::edge_reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1166 
1167 template LIBMESH_EXPORT void FE<2,RAVIART_THOMAS>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1168 template LIBMESH_EXPORT void FE<2,RAVIART_THOMAS>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1169 template LIBMESH_EXPORT void FE<2,RAVIART_THOMAS>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1170 template LIBMESH_EXPORT void FE<2,RAVIART_THOMAS>::edge_reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1171 
1172 template LIBMESH_EXPORT void FE<2,L2_RAVIART_THOMAS>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1173 template LIBMESH_EXPORT void FE<2,L2_RAVIART_THOMAS>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1174 template LIBMESH_EXPORT void FE<2,L2_RAVIART_THOMAS>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1175 template LIBMESH_EXPORT void FE<2,L2_RAVIART_THOMAS>::edge_reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1176 
1177 template LIBMESH_EXPORT void FE<3,NEDELEC_ONE>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1178 template LIBMESH_EXPORT void FE<3,NEDELEC_ONE>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1179 template LIBMESH_EXPORT void FE<3,NEDELEC_ONE>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1180 template LIBMESH_EXPORT void FE<3,NEDELEC_ONE>::edge_reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1181 
1182 template LIBMESH_EXPORT void FE<3,RAVIART_THOMAS>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1183 template LIBMESH_EXPORT void FE<3,RAVIART_THOMAS>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1184 template LIBMESH_EXPORT void FE<3,RAVIART_THOMAS>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1185 template LIBMESH_EXPORT void FE<3,RAVIART_THOMAS>::edge_reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1186 
1187 template LIBMESH_EXPORT void FE<3,L2_RAVIART_THOMAS>::reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1188 template LIBMESH_EXPORT void FE<3,L2_RAVIART_THOMAS>::side_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1189 template LIBMESH_EXPORT void FE<3,L2_RAVIART_THOMAS>::edge_map(Elem const *, Elem const *, const unsigned int, const std::vector<Point> &, std::vector<Point> &);
1190 template LIBMESH_EXPORT void FE<3,L2_RAVIART_THOMAS>::edge_reinit(Elem const *, unsigned int, Real, const std::vector<Point> * const, const std::vector<Real> * const);
1191 
1192 // Intel 9.1 complained it needed this in devel mode.
1193 //template LIBMESH_EXPORT void FE<2,XYZ>::init_face_shape_functions(const std::vector<Point> &, const Elem *);
1194 //template LIBMESH_EXPORT void FE<3,XYZ>::init_face_shape_functions(const std::vector<Point> &, const Elem *);
1195 
1196 } // namespace libMesh
Real(* shape_deriv_ptr)(const FEType fet, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
Typedef for pointer to a function that returns FE shape function derivative values.
Definition: fe_interface.h:623
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
void init_edge_shape_functions(const std::vector< Point > &qp, const Elem *edge)
Same as before, but for an edge.
Definition: fe_boundary.C:565
Real dzdeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected...
Definition: fe_map.h:711
std::vector< std::vector< Real > > dpsideta_map
Map for the derivative of the side function, d(psi)/d(eta).
Definition: fe_map.h:945
std::vector< std::vector< Real > > psi_map
Map for the side shape functions, psi.
Definition: fe_map.h:933
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
std::vector< std::vector< Real > > d2psideta2_map
Map for the second derivatives (in eta) of the side shape functions.
Definition: fe_map.h:968
virtual void side_map(const Elem *elem, const Elem *side, const unsigned int s, const std::vector< Point > &reference_side_points, std::vector< Point > &reference_points) override
Computes the reference space quadrature points on the side of an element based on the side quadrature...
Definition: fe_boundary.C:336
static shape_ptr shape_function(const unsigned int dim, const FEType &fe_t, const ElemType t)
const Elem * interior_parent() const
Definition: elem.C:994
bool calculate_dxyz
Should we calculate mapping gradients?
Definition: fe_map.h:1016
void compute_edge_map(int dim, const std::vector< Real > &qw, const Elem *side)
Same as before, but for an edge.
Definition: fe_boundary.C:1005
unsigned int dim
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true)
Definition: fe_map.C:1626
Real dxdeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected...
Definition: fe_map.h:695
Real dzdxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected...
Definition: fe_map.h:687
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
unsigned int p_level() const
Definition: elem.h:2945
ADRealEigenVector< T, D, asd > sqrt(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:53
The libMesh namespace provides an interface to certain functionality in the library.
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy=false)=0
const Number zero
.
Definition: libmesh.h:280
std::vector< std::vector< Real > > d2psidxideta_map
Map for the second (cross) derivatives in xi, eta of the side shape functions.
Definition: fe_map.h:961
std::vector< RealGradient > d2xyzdxideta_map
Vector of mixed second partial derivatives in xi-eta: d^2(x)/d(xi)d(eta) d^2(y)/d(xi)d(eta) d^2(z)/d(...
Definition: fe_map.h:772
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
This is at the core of this class.
Definition: fe.C:142
std::vector< RealGradient > dxyzdxi_map
Vector of partial derivatives: d(x)/d(xi), d(y)/d(xi), d(z)/d(xi)
Definition: fe_map.h:746
TypeVector< T > unit() const
Definition: type_vector.h:1120
std::vector< RealGradient > d2xyzdeta2_map
Vector of second partial derivatives in eta: d^2(x)/d(eta)^2.
Definition: fe_map.h:778
REINIT_ERROR(1, RAVIART_THOMAS, reinit)
Definition: fe_boundary.C:101
SIDEMAP_ERROR(1, NEDELEC_ONE, side_map)
Definition: fe_boundary.C:100
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:515
std::vector< RealGradient > d2xyzdxi2_map
Vector of second partial derivatives in xi: d^2(x)/d(xi)^2, d^2(y)/d(xi)^2, d^2(z)/d(xi)^2.
Definition: fe_map.h:766
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const =0
LIBMESH_ERRORS_IN_LOW_D(CLOUGH)
Definition: fe_boundary.C:73
bool calculate_d2xyz
Should we calculate mapping hessians?
Definition: fe_map.h:1023
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const =0
Similar to Elem::local_side_node(), but instead of a side id, takes an edge id and a node id on that ...
libmesh_assert(ctx)
bool calculate_xyz
Should we calculate physical point locations?
Definition: fe_map.h:1011
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
Definition: type_vector.h:906
std::vector< Real > curvatures
The mean curvature (= one half the sum of the principal curvatures) on the boundary at the quadrature...
Definition: fe_map.h:988
std::vector< std::vector< Real > > d2psidxi2_map
Map for the second derivatives (in xi) of the side shape functions.
Definition: fe_map.h:954
static shape_deriv_ptr shape_deriv_function(const unsigned int dim, const FEType &fe_t, const ElemType t)
std::vector< Real > JxW
Jacobian*Weight values at quadrature points.
Definition: fe_map.h:1000
std::vector< Point > normals
Normal vectors on boundary at quadrature points.
Definition: fe_map.h:980
std::vector< Point > xyz
The spatial locations of the quadrature points.
Definition: fe_map.h:740
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2407
Real(* shape_ptr)(const FEType fe_t, const Elem *elem, const unsigned int i, const Point &p, const bool add_p_level)
Typedef for pointer to a function that returns FE shape function values.
Definition: fe_interface.h:516
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real dydeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected...
Definition: fe_map.h:703
static shape_second_deriv_ptr shape_second_deriv_function(const unsigned int dim, const FEType &fe_t, const ElemType t)
virtual void edge_map(const Elem *elem, const Elem *edge, const unsigned int e, const std::vector< Point > &reference_edge_points, std::vector< Point > &reference_points)
Computes the reference space quadrature points on the side of an element based on the edge quadrature...
Definition: fe_boundary.C:387
std::vector< std::vector< Real > > dpsidxi_map
Map for the derivative of the side functions, d(psi)/d(xi).
Definition: fe_map.h:939
std::vector< std::vector< Point > > tangents
Tangent vectors on boundary at quadrature points.
Definition: fe_map.h:975
std::vector< RealGradient > dxyzdeta_map
Vector of partial derivatives: d(x)/d(eta), d(y)/d(eta), d(z)/d(eta)
Definition: fe_map.h:752
static Point map_deriv(const unsigned int dim, const Elem *elem, const unsigned int j, const Point &reference_point)
Definition: fe_map.C:2101
Real(* shape_second_deriv_ptr)(const FEType fet, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
Typedef for pointer to a function that returns FE shape function second derivative values...
Definition: fe_interface.h:727
Real dxdxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected...
Definition: fe_map.h:671
void determine_calculations()
Determine which values are to be calculated.
Definition: fe_map.h:648
Real dydxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected...
Definition: fe_map.h:679
virtual void compute_face_map(int dim, const std::vector< Real > &qw, const Elem *side)
Same as compute_map, but for a side.
Definition: fe_boundary.C:641
FEFamily
defines an enum for finite element families.
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i)=0
virtual Order default_order() const =0
virtual void edge_reinit(const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
Reinitializes all the physical element-dependent data based on the edge.
Definition: fe_boundary.C:237
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:2299
const Point & point(const unsigned int i) const
Definition: elem.h:2277
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
void init_face_shape_functions(const std::vector< Point > &qp, const Elem *side)
Initializes the reference to physical element map for a side.
Definition: fe_boundary.C:437
static FEFamily map_fe_type(const Elem &elem)
Definition: fe_map.C:45