libMesh::FEInterface Class Reference

Interface class which provides access to FE functions. More...

#include <fe_interface.h>

Public Types

typedef unsigned int(* n_dofs_at_node_ptr) (const ElemType, const Order, const unsigned int)
 

Public Member Functions

virtual ~FEInterface ()
 
template<>
void shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p, Real &phi)
 
template<>
void shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, Real &phi)
 
template<>
void shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p, RealGradient &phi)
 
template<>
void shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, RealGradient &phi)
 

Static Public Member Functions

static unsigned int n_shape_functions (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
static unsigned int n_dofs (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
static unsigned int n_dofs_at_node (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
 
static n_dofs_at_node_ptr n_dofs_at_node_function (const unsigned int dim, const FEType &fe_t)
 
static unsigned int n_dofs_per_elem (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
static void dofs_on_side (const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int s, std::vector< unsigned int > &di)
 
static void dofs_on_edge (const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int e, std::vector< unsigned int > &di)
 
static void nodal_soln (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
 
static Point map (unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p)
 
static Point inverse_map (const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
 
static void inverse_map (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &physical_points, std::vector< Point > &reference_points, const Real tolerance=TOLERANCE, const bool secure=true)
 
static bool on_reference_element (const Point &p, const ElemType t, const Real eps=TOLERANCE)
 
static Real shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
 
static Real shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p)
 
template<typename OutputType >
static void shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p, OutputType &phi)
 
template<typename OutputType >
static void shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p, OutputType &phi)
 
static void compute_data (const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
 
static void compute_constraints (DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
 
static void compute_periodic_constraints (DofConstraints &constraints, DofMap &dof_map, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const unsigned int variable_number, const Elem *elem)
 
static unsigned int max_order (const FEType &fe_t, const ElemType &el_t)
 
static bool extra_hanging_dofs (const FEType &fe_t)
 
static FEFieldType field_type (const FEType &fe_type)
 
static FEFieldType field_type (const FEFamily &fe_family)
 
static unsigned int n_vec_dim (const MeshBase &mesh, const FEType &fe_type)
 

Private Member Functions

 FEInterface ()
 

Static Private Member Functions

static bool is_InfFE_elem (const ElemType et)
 
static unsigned int ifem_n_shape_functions (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
static unsigned int ifem_n_dofs (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
static unsigned int ifem_n_dofs_at_node (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
 
static unsigned int ifem_n_dofs_per_elem (const unsigned int dim, const FEType &fe_t, const ElemType t)
 
static void ifem_nodal_soln (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
 
static Point ifem_map (const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p)
 
static Point ifem_inverse_map (const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
 
static void ifem_inverse_map (const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Point > &physical_points, std::vector< Point > &reference_points, const Real tolerance=TOLERANCE, const bool secure=true)
 
static bool ifem_on_reference_element (const Point &p, const ElemType t, const Real eps)
 
static Real ifem_shape (const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
 
static Real ifem_shape (const unsigned int dim, const FEType &fe_t, const Elem *elem, const unsigned int i, const Point &p)
 
static void ifem_compute_data (const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
 

Detailed Description

Interface class which provides access to FE functions.

This class provides an encapsulated access to all static public member functions of finite element classes. Using this class, one need not worry about the correct finite element class.

Author
Daniel Dreyer
Date
2002-2007

Definition at line 64 of file fe_interface.h.

Member Typedef Documentation

◆ n_dofs_at_node_ptr

typedef unsigned int(* libMesh::FEInterface::n_dofs_at_node_ptr) (const ElemType, const Order, const unsigned int)

Definition at line 114 of file fe_interface.h.

Constructor & Destructor Documentation

◆ FEInterface()

libMesh::FEInterface::FEInterface ( )
private

Empty constructor. Do not create an object of this type.

Definition at line 35 of file fe_interface.C.

36 {
37  libmesh_error_msg("ERROR: Do not define an object of this type.");
38 }

◆ ~FEInterface()

virtual libMesh::FEInterface::~FEInterface ( )
inlinevirtual

Destructor.

Definition at line 78 of file fe_interface.h.

78 {}

Member Function Documentation

◆ compute_constraints()

void libMesh::FEInterface::compute_constraints ( DofConstraints constraints,
DofMap dof_map,
const unsigned int  variable_number,
const Elem elem 
)
static

Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to variable number var_number.

Definition at line 873 of file fe_interface.C.

References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::FE< Dim, T >::compute_constraints(), libMesh::Elem::dim(), libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::L2_HIERARCHIC, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::SZABAB, and libMesh::DofMap::variable_type().

877 {
878  libmesh_assert(elem);
879 
880  const FEType & fe_t = dof_map.variable_type(variable_number);
881 
882  switch (elem->dim())
883  {
884  case 0:
885  case 1:
886  {
887  // No constraints in 0D/1D.
888  return;
889  }
890 
891 
892  case 2:
893  {
894  switch (fe_t.family)
895  {
896  case CLOUGH:
898  dof_map,
899  variable_number,
900  elem); return;
901 
902  case HERMITE:
904  dof_map,
905  variable_number,
906  elem); return;
907 
908  case LAGRANGE:
910  dof_map,
911  variable_number,
912  elem); return;
913 
914  case HIERARCHIC:
916  dof_map,
917  variable_number,
918  elem); return;
919 
920  case L2_HIERARCHIC:
922  dof_map,
923  variable_number,
924  elem); return;
925 
926  case LAGRANGE_VEC:
928  dof_map,
929  variable_number,
930  elem); return;
931 
932 
933 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
934  case SZABAB:
936  dof_map,
937  variable_number,
938  elem); return;
939 
940  case BERNSTEIN:
942  dof_map,
943  variable_number,
944  elem); return;
945 
946 #endif
947  default:
948  return;
949  }
950  }
951 
952 
953  case 3:
954  {
955  switch (fe_t.family)
956  {
957  case HERMITE:
959  dof_map,
960  variable_number,
961  elem); return;
962 
963  case LAGRANGE:
965  dof_map,
966  variable_number,
967  elem); return;
968 
969  case HIERARCHIC:
971  dof_map,
972  variable_number,
973  elem); return;
974 
975  case L2_HIERARCHIC:
977  dof_map,
978  variable_number,
979  elem); return;
980 
981  case LAGRANGE_VEC:
983  dof_map,
984  variable_number,
985  elem); return;
986 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
987  case SZABAB:
989  dof_map,
990  variable_number,
991  elem); return;
992 
993  case BERNSTEIN:
995  dof_map,
996  variable_number,
997  elem); return;
998 
999 #endif
1000  default:
1001  return;
1002  }
1003  }
1004 
1005 
1006  default:
1007  libmesh_error_msg("Invalid dimension = " << elem->dim());
1008  }
1009 }
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)

◆ compute_data()

void libMesh::FEInterface::compute_data ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
FEComputeData data 
)
static

Lets the appropriate child of FEBase compute the requested data for the input specified in data, and sets the values in data. See this as a generalization of shape(). With infinite elements disabled, computes values for all shape functions of elem evaluated at p.

Note
On a p-refined element, fe_t.order should be the base order of the element.

Definition at line 837 of file fe_interface.C.

References data, ifem_compute_data(), is_InfFE_elem(), n_dofs(), libMesh::FEType::order, libMesh::Elem::p_level(), shape(), and libMesh::Elem::type().

Referenced by libMesh::MeshFunction::discontinuous_value(), libMesh::MeshFunction::operator()(), and libMesh::System::point_value().

841 {
842 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
843 
844  if (elem && is_InfFE_elem(elem->type()))
845  {
846  data.init();
847  ifem_compute_data(dim, fe_t, elem, data);
848  return;
849  }
850 
851 #endif
852 
853  FEType p_refined = fe_t;
854  p_refined.order = static_cast<Order>(p_refined.order + elem->p_level());
855 
856  const unsigned int n_dof = n_dofs (dim, p_refined, elem->type());
857  const Point & p = data.p;
858  data.shape.resize(n_dof);
859 
860  // set default values for all the output fields
861  data.init();
862 
863  for (unsigned int n=0; n<n_dof; n++)
864  data.shape[n] = shape(dim, p_refined, elem, n, p);
865 
866  return;
867 }
static unsigned int n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:454
static void ifem_compute_data(const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
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:657
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45
IterBase * data

◆ compute_periodic_constraints()

void libMesh::FEInterface::compute_periodic_constraints ( DofConstraints constraints,
DofMap dof_map,
const PeriodicBoundaries boundaries,
const MeshBase mesh,
const PointLocatorBase point_locator,
const unsigned int  variable_number,
const Elem elem 
)
static

Computes the constraint matrix contributions (for periodic boundary conditions) corresponding to variable number var_number.

Definition at line 1017 of file fe_interface.C.

References libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), and mesh.

1024 {
1025  // No element-specific optimizations currently exist
1027  dof_map,
1028  boundaries,
1029  mesh,
1030  point_locator,
1031  variable_number,
1032  elem);
1033 }
MeshBase & mesh
static void compute_periodic_constraints(DofConstraints &constraints, DofMap &dof_map, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const unsigned int variable_number, const Elem *elem)
Definition: fe_base.C:1650

◆ dofs_on_edge()

void libMesh::FEInterface::dofs_on_edge ( const Elem *const  elem,
const unsigned int  dim,
const FEType fe_t,
unsigned int  e,
std::vector< unsigned int > &  di 
)
static

Fills the vector di with the local degree of freedom indices associated with edge e of element elem Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the base order of the element.

Definition at line 536 of file fe_interface.C.

References libMesh::FEType::order.

Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), and libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()().

541 {
542  const Order o = fe_t.order;
543 
544  void_fe_with_vec_switch(dofs_on_edge(elem, o, e, di));
545 }
static void dofs_on_edge(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int e, std::vector< unsigned int > &di)
Definition: fe_interface.C:536

◆ dofs_on_side()

void libMesh::FEInterface::dofs_on_side ( const Elem *const  elem,
const unsigned int  dim,
const FEType fe_t,
unsigned int  s,
std::vector< unsigned int > &  di 
)
static

Fills the vector di with the local degree of freedom indices associated with side s of element elem Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the base order of the element.

Definition at line 523 of file fe_interface.C.

References libMesh::FEType::order.

Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_periodic_constraints(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_proj_constraints(), and libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()().

528 {
529  const Order o = fe_t.order;
530 
531  void_fe_with_vec_switch(dofs_on_side(elem, o, s, di));
532 }
static void dofs_on_side(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int s, std::vector< unsigned int > &di)
Definition: fe_interface.C:523

◆ extra_hanging_dofs()

bool libMesh::FEInterface::extra_hanging_dofs ( const FEType fe_t)
static
Returns
true if separate degrees of freedom must be allocated for vertex DoFs and edge/face DoFs at a hanging node.

Definition at line 1412 of file fe_interface.C.

References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::L2_HIERARCHIC, libMesh::L2_LAGRANGE, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::NEDELEC_ONE, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.

Referenced by libMesh::DofMap::_dof_indices(), libMesh::DofMap::old_dof_indices(), and libMesh::DofMap::reinit().

1413 {
1414  switch (fe_t.family)
1415  {
1416  case LAGRANGE:
1417  case L2_LAGRANGE:
1418  case MONOMIAL:
1419  case L2_HIERARCHIC:
1420  case XYZ:
1421  case SUBDIVISION:
1422  case LAGRANGE_VEC:
1423  case NEDELEC_ONE:
1424  return false;
1425  case CLOUGH:
1426  case HERMITE:
1427  case HIERARCHIC:
1428 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
1429  case BERNSTEIN:
1430  case SZABAB:
1431 #endif
1432  default:
1433  return true;
1434  }
1435 }

◆ field_type() [1/2]

FEFieldType libMesh::FEInterface::field_type ( const FEType fe_type)
static

◆ field_type() [2/2]

FEFieldType libMesh::FEInterface::field_type ( const FEFamily fe_family)
static
Returns
The number of components of a vector-valued element. Scalar-valued elements return 1.

Definition at line 1442 of file fe_interface.C.

References libMesh::LAGRANGE_VEC, libMesh::NEDELEC_ONE, libMesh::TYPE_SCALAR, and libMesh::TYPE_VECTOR.

1443 {
1444  switch (fe_family)
1445  {
1446  case LAGRANGE_VEC:
1447  case NEDELEC_ONE:
1448  return TYPE_VECTOR;
1449  default:
1450  return TYPE_SCALAR;
1451  }
1452 }

◆ ifem_compute_data()

void libMesh::FEInterface::ifem_compute_data ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
FEComputeData data 
)
staticprivate

Definition at line 808 of file fe_interface_inf_fe.C.

References libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), data, libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::LAGRANGE, libMesh::LEGENDRE, and libMesh::FEType::radial_family.

Referenced by compute_data().

812 {
813  switch (dim)
814  {
815  // 1D
816  case 1:
817  {
818  switch (fe_t.radial_family)
819  {
820  /*
821  * For no derivatives (and local coordinates, as
822  * given in \p p) the infinite element shapes
823  * are independent of mapping type
824  */
825  case INFINITE_MAP:
827  break;
828 
829  case JACOBI_20_00:
831  break;
832 
833  case JACOBI_30_00:
835  break;
836 
837  case LEGENDRE:
839  break;
840 
841  case LAGRANGE:
843  break;
844 
845  default:
846  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
847  }
848 
849  break;
850  }
851 
852 
853  // 2D
854  case 2:
855  {
856  switch (fe_t.radial_family)
857  {
858  case INFINITE_MAP:
860  break;
861 
862  case JACOBI_20_00:
864  break;
865 
866  case JACOBI_30_00:
868  break;
869 
870  case LEGENDRE:
872  break;
873 
874  case LAGRANGE:
876  break;
877 
878  default:
879  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
880  }
881 
882  break;
883  }
884 
885 
886  // 3D
887  case 3:
888  {
889  switch (fe_t.radial_family)
890  {
891  case INFINITE_MAP:
893  break;
894 
895  case JACOBI_20_00:
897  break;
898 
899  case JACOBI_30_00:
901  break;
902 
903  case LEGENDRE:
905  break;
906 
907  case LAGRANGE:
909  break;
910 
911  default:
912  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
913  }
914 
915  break;
916  }
917 
918 
919  default:
920  libmesh_error_msg("Invalid dim = " << dim);
921  break;
922  }
923 }
static void compute_data(const FEType &fe_t, const Elem *inf_elem, FEComputeData &data)
IterBase * data

◆ ifem_inverse_map() [1/2]

Point libMesh::FEInterface::ifem_inverse_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const Point p,
const Real  tolerance = TOLERANCE,
const bool  secure = true 
)
staticprivate

Definition at line 447 of file fe_interface_inf_fe.C.

References libMesh::CARTESIAN, libMesh::ELLIPSOIDAL, libMesh::FEType::inf_map, libMesh::InfFE< Dim, T_radial, T_map >::inverse_map(), and libMesh::SPHERICAL.

Referenced by inverse_map().

453 {
454  switch (dim)
455  {
456  // 1D
457  case 1:
458  {
459  switch (fe_t.inf_map)
460  {
461  case CARTESIAN:
462  return InfFE<1,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure);
463 
464  case SPHERICAL:
465  case ELLIPSOIDAL:
466  libmesh_not_implemented_msg("ERROR: Spherical and Ellipsoidal IFEMs not (yet) implemented.");
467 
468  /*
469  case SPHERICAL:
470  return InfFE<1,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance);
471 
472  case ELLIPSOIDAL:
473  return InfFE<1,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance);
474  */
475 
476  default:
477  libmesh_error_msg("Invalid map = " << fe_t.inf_map);
478  }
479  }
480 
481 
482  // 2D
483  case 2:
484  {
485  switch (fe_t.inf_map)
486  {
487  case CARTESIAN:
488  return InfFE<2,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure);
489 
490  case SPHERICAL:
491  case ELLIPSOIDAL:
492  libmesh_not_implemented_msg("ERROR: Spherical and Ellipsoidal IFEMs not (yet) implemented.");
493 
494  /*
495  case SPHERICAL:
496  return InfFE<2,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance);
497 
498  case ELLIPSOIDAL:
499  return InfFE<2,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance);
500  */
501 
502  default:
503  libmesh_error_msg("Invalid map = " << fe_t.inf_map);
504  }
505  }
506 
507 
508  // 3D
509  case 3:
510  {
511  switch (fe_t.inf_map)
512  {
513  case CARTESIAN:
514  return InfFE<3,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure);
515 
516  case SPHERICAL:
517  case ELLIPSOIDAL:
518  libmesh_not_implemented_msg("ERROR: Spherical and Ellipsoidal IFEMs not (yet) implemented.");
519 
520  /*
521  case SPHERICAL:
522  return InfFE<3,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance);
523 
524  case ELLIPSOIDAL:
525  return InfFE<3,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance);
526  */
527 
528  default:
529  libmesh_error_msg("Invalid map = " << fe_t.inf_map);
530  }
531  }
532 
533  default:
534  libmesh_error_msg("Invalid dim = " << dim);
535  }
536 }
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: inf_fe_map.C:90

◆ ifem_inverse_map() [2/2]

void libMesh::FEInterface::ifem_inverse_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const std::vector< Point > &  physical_points,
std::vector< Point > &  reference_points,
const Real  tolerance = TOLERANCE,
const bool  secure = true 
)
staticprivate

Definition at line 540 of file fe_interface_inf_fe.C.

References libMesh::CARTESIAN, libMesh::FEType::inf_map, and libMesh::InfFE< Dim, T_radial, T_map >::inverse_map().

547 {
548  switch (dim)
549  {
550  // 1D
551  case 1:
552  {
553  switch (fe_t.inf_map)
554  {
555  case CARTESIAN:
556  InfFE<1,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure);
557  return;
558 
559  default:
560  libmesh_error_msg("Invalid map = " << fe_t.inf_map);
561  }
562  }
563 
564 
565  // 2D
566  case 2:
567  {
568  switch (fe_t.inf_map)
569  {
570  case CARTESIAN:
571  InfFE<2,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure);
572  return;
573 
574  default:
575  libmesh_error_msg("Invalid map = " << fe_t.inf_map);
576  }
577  }
578 
579 
580  // 3D
581  case 3:
582  {
583  switch (fe_t.inf_map)
584  {
585  case CARTESIAN:
586  InfFE<3,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure);
587  return;
588 
589  default:
590  libmesh_error_msg("Invalid map = " << fe_t.inf_map);
591  }
592  }
593 
594  default:
595  libmesh_error_msg("Invalid dim = " << dim);
596  }
597 }
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: inf_fe_map.C:90

◆ ifem_map()

Point libMesh::FEInterface::ifem_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const Point p 
)
staticprivate

Definition at line 416 of file fe_interface_inf_fe.C.

References libMesh::CARTESIAN, libMesh::ELLIPSOIDAL, libMesh::FEType::inf_map, libMesh::InfFE< Dim, T_radial, T_map >::map(), and libMesh::SPHERICAL.

Referenced by map().

420 {
421  switch (fe_t.inf_map)
422  {
423  case CARTESIAN:
424  {
425  switch (dim)
426  {
427  case 1:
428  return InfFE<1,JACOBI_20_00,CARTESIAN>::map(elem, p);
429  case 2:
430  return InfFE<2,JACOBI_20_00,CARTESIAN>::map(elem, p);
431  case 3:
432  return InfFE<3,JACOBI_20_00,CARTESIAN>::map(elem, p);
433  default:
434  libmesh_error_msg("Invalid dim = " << dim);
435  }
436  }
437  case SPHERICAL:
438  case ELLIPSOIDAL:
439  libmesh_not_implemented_msg("ERROR: Spherical and Ellipsoidal IFEMs not (yet) implemented.");
440  default:
441  libmesh_error_msg("Invalid map = " << fe_t.inf_map);
442  }
443 }
static Point map(const Elem *inf_elem, const Point &reference_point)
Definition: inf_fe_map.C:40

◆ ifem_n_dofs()

unsigned int libMesh::FEInterface::ifem_n_dofs ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
staticprivate

Definition at line 71 of file fe_interface_inf_fe.C.

References libMesh::InfFE< Dim, T_radial, T_map >::n_dofs().

Referenced by n_dofs().

74 {
75  switch (dim)
76  {
77  // 1D
78  case 1:
79  /*
80  * Since InfFE<Dim,T_radial,T_map>::n_dofs(...)
81  * is actually independent of T_radial and T_map, we can use
82  * just any T_radial and T_map
83  */
85 
86  // 2D
87  case 2:
89 
90  // 3D
91  case 3:
93 
94  default:
95  libmesh_error_msg("Unsupported dim = " << dim);
96  }
97 }
static unsigned int n_dofs(const FEType &fet, const ElemType inf_elem_type)
Definition: inf_fe_static.C:55

◆ ifem_n_dofs_at_node()

unsigned int libMesh::FEInterface::ifem_n_dofs_at_node ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  n 
)
staticprivate

Definition at line 102 of file fe_interface_inf_fe.C.

References libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_at_node().

Referenced by n_dofs_at_node().

106 {
107  switch (dim)
108  {
109  // 1D
110  case 1:
111  /*
112  * Since InfFE<Dim,T_radial,T_map>::n_dofs_at_node(...)
113  * is actually independent of T_radial and T_map, we can use
114  * just any T_radial and T_map
115  */
117 
118  // 2D
119  case 2:
121 
122  // 3D
123  case 3:
125 
126  default:
127  libmesh_error_msg("Unsupported dim = " << dim);
128  }
129 }
static unsigned int n_dofs_at_node(const FEType &fet, const ElemType inf_elem_type, const unsigned int n)
Definition: inf_fe_static.C:72

◆ ifem_n_dofs_per_elem()

unsigned int libMesh::FEInterface::ifem_n_dofs_per_elem ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
staticprivate

Definition at line 135 of file fe_interface_inf_fe.C.

References libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem().

Referenced by n_dofs_per_elem().

138 {
139  switch (dim)
140  {
141  // 1D
142  case 1:
143  /*
144  * Since InfFE<Dim,T_radial,T_map>::n_dofs(...)
145  * is actually independent of T_radial and T_map, we can use
146  * just any T_radial and T_map
147  */
149 
150  // 2D
151  case 2:
153 
154  // 3D
155  case 3:
157 
158  default:
159  libmesh_error_msg("Unsupported dim = " << dim);
160  }
161 }
static unsigned int n_dofs_per_elem(const FEType &fet, const ElemType inf_elem_type)

◆ ifem_n_shape_functions()

unsigned int libMesh::FEInterface::ifem_n_shape_functions ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
staticprivate

Definition at line 39 of file fe_interface_inf_fe.C.

References libMesh::InfFE< Dim, T_radial, T_map >::n_shape_functions().

Referenced by n_shape_functions().

42 {
43  switch (dim)
44  {
45  // 1D
46  case 1:
47  /*
48  * Since InfFE<Dim,T_radial,T_map>::n_shape_functions(...)
49  * is actually independent of T_radial and T_map, we can use
50  * just any T_radial and T_map
51  */
53 
54  // 2D
55  case 2:
57 
58  // 3D
59  case 3:
61 
62  default:
63  libmesh_error_msg("Unsupported dim = " << dim);
64  }
65 }
virtual unsigned int n_shape_functions() const override
Definition: inf_fe.h:460

◆ ifem_nodal_soln()

void libMesh::FEInterface::ifem_nodal_soln ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const std::vector< Number > &  elem_soln,
std::vector< Number > &  nodal_soln 
)
staticprivate

Definition at line 166 of file fe_interface_inf_fe.C.

References libMesh::CARTESIAN, libMesh::FEType::inf_map, libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::LAGRANGE, libMesh::LEGENDRE, nodal_soln(), libMesh::InfFE< Dim, T_radial, T_map >::nodal_soln(), and libMesh::FEType::radial_family.

Referenced by nodal_soln().

171 {
172  switch (dim)
173  {
174 
175  // 1D
176  case 1:
177  {
178  switch (fe_t.radial_family)
179  {
180  case INFINITE_MAP:
181  libmesh_error_msg("ERROR: INFINITE_MAP is not a valid shape family for radial approximation.");
182 
183  case JACOBI_20_00:
184  {
185  switch (fe_t.inf_map)
186  {
187  case CARTESIAN:
188  {
190  break;
191  }
192  default:
193  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
194  }
195  break;
196  }
197 
198  case JACOBI_30_00:
199  {
200  switch (fe_t.inf_map)
201  {
202  case CARTESIAN:
203  {
205  break;
206  }
207  default:
208  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
209  }
210  break;
211  }
212 
213  case LEGENDRE:
214  {
215  switch (fe_t.inf_map)
216  {
217  case CARTESIAN:
218  {
219  InfFE<1,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
220  break;
221  }
222  default:
223  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
224  }
225  break;
226  }
227 
228  case LAGRANGE:
229  {
230  switch (fe_t.inf_map)
231  {
232  case CARTESIAN:
233  {
234  InfFE<1,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
235  break;
236  }
237  default:
238  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
239  }
240  break;
241  }
242 
243  default:
244  libmesh_error_msg("ERROR: Bad FEType.radial_family= " << fe_t.radial_family);
245  }
246 
247  break;
248  }
249 
250 
251 
252 
253  // 2D
254  case 2:
255  {
256  switch (fe_t.radial_family)
257  {
258  case INFINITE_MAP:
259  libmesh_error_msg("ERROR: INFINITE_MAP is not a valid shape family for radial approximation.");
260 
261  case JACOBI_20_00:
262  {
263  switch (fe_t.inf_map)
264  {
265  case CARTESIAN:
266  {
268  break;
269  }
270  default:
271  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
272  }
273  break;
274  }
275 
276  case JACOBI_30_00:
277  {
278  switch (fe_t.inf_map)
279  {
280  case CARTESIAN:
281  {
283  break;
284  }
285  default:
286  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
287  }
288  break;
289  }
290 
291  case LEGENDRE:
292  {
293  switch (fe_t.inf_map)
294  {
295  case CARTESIAN:
296  {
297  InfFE<2,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
298  break;
299  }
300  default:
301  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
302  }
303  break;
304  }
305 
306  case LAGRANGE:
307  {
308  switch (fe_t.inf_map)
309  {
310  case CARTESIAN:
311  {
312  InfFE<2,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
313  break;
314  }
315  default:
316  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
317  }
318  break;
319  }
320 
321  default:
322  libmesh_error_msg("ERROR: Bad FEType.radial_family= " << fe_t.radial_family);
323  }
324 
325  break;
326  }
327 
328 
329 
330 
331  // 3D
332  case 3:
333  {
334  switch (fe_t.radial_family)
335  {
336  case INFINITE_MAP:
337  libmesh_error_msg("ERROR: INFINITE_MAP is not a valid shape family for radial approximation.");
338 
339  case JACOBI_20_00:
340  {
341  switch (fe_t.inf_map)
342  {
343  case CARTESIAN:
344  {
346  break;
347  }
348  default:
349  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
350  }
351  break;
352  }
353 
354  case JACOBI_30_00:
355  {
356  switch (fe_t.inf_map)
357  {
358  case CARTESIAN:
359  {
361  break;
362  }
363  default:
364  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
365  }
366  break;
367  }
368 
369  case LEGENDRE:
370  {
371  switch (fe_t.inf_map)
372  {
373  case CARTESIAN:
374  {
375  InfFE<3,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
376  break;
377  }
378  default:
379  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
380  }
381  break;
382  }
383 
384  case LAGRANGE:
385  {
386  switch (fe_t.inf_map)
387  {
388  case CARTESIAN:
389  {
390  InfFE<3,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln);
391  break;
392  }
393  default:
394  libmesh_error_msg("ERROR: Spherical & Ellipsoidal IFEMs not implemented.");
395  }
396  break;
397  }
398 
399 
400 
401  default:
402  libmesh_error_msg("ERROR: Bad FEType.radial_family= " << fe_t.radial_family);
403  }
404 
405  break;
406  }
407 
408  default:
409  libmesh_error_msg("Invalid dim = " << dim);
410  }
411 }
static void nodal_soln(const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
Definition: fe_interface.C:550
static void nodal_soln(const FEType &fet, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)

◆ ifem_on_reference_element()

bool libMesh::FEInterface::ifem_on_reference_element ( const Point p,
const ElemType  t,
const Real  eps 
)
staticprivate

Definition at line 602 of file fe_interface_inf_fe.C.

References libMesh::FEAbstract::on_reference_element().

605 {
606  return FEBase::on_reference_element(p,t,eps);
607 }
static bool on_reference_element(const Point &p, const ElemType t, const Real eps=TOLERANCE)
Definition: fe_abstract.C:581

◆ ifem_shape() [1/2]

Real libMesh::FEInterface::ifem_shape ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  i,
const Point p 
)
staticprivate

Definition at line 612 of file fe_interface_inf_fe.C.

References libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::LAGRANGE, libMesh::LEGENDRE, libMesh::FEType::radial_family, and libMesh::InfFE< Dim, T_radial, T_map >::shape().

Referenced by shape().

617 {
618  switch (dim)
619  {
620  // 1D
621  case 1:
622  {
623  switch (fe_t.radial_family)
624  {
625  /*
626  * For no derivatives (and local coordinates, as
627  * given in \p p) the infinite element shapes
628  * are independent of mapping type
629  */
630  case INFINITE_MAP:
631  return InfFE<1,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p);
632 
633  case JACOBI_20_00:
634  return InfFE<1,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p);
635 
636  case JACOBI_30_00:
637  return InfFE<1,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p);
638 
639  case LEGENDRE:
640  return InfFE<1,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p);
641 
642  case LAGRANGE:
643  return InfFE<1,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p);
644 
645  default:
646  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
647  }
648  }
649 
650 
651  // 2D
652  case 2:
653  {
654  switch (fe_t.radial_family)
655  {
656  case INFINITE_MAP:
657  return InfFE<2,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p);
658 
659  case JACOBI_20_00:
660  return InfFE<2,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p);
661 
662  case JACOBI_30_00:
663  return InfFE<2,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p);
664 
665  case LEGENDRE:
666  return InfFE<2,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p);
667 
668  case LAGRANGE:
669  return InfFE<2,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p);
670 
671  default:
672  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
673  }
674  }
675 
676 
677  // 3D
678  case 3:
679  {
680  switch (fe_t.radial_family)
681  {
682  case INFINITE_MAP:
683  return InfFE<3,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p);
684 
685  case JACOBI_20_00:
686  return InfFE<3,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p);
687 
688  case JACOBI_30_00:
689  return InfFE<3,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p);
690 
691  case LEGENDRE:
692  return InfFE<3,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p);
693 
694  case LAGRANGE:
695  return InfFE<3,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p);
696 
697  default:
698  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
699  }
700  }
701 
702  default:
703  libmesh_error_msg("Invalid dim = " << dim);
704  }
705 }
static Real shape(const FEType &fet, const ElemType t, const unsigned int i, const Point &p)

◆ ifem_shape() [2/2]

Real libMesh::FEInterface::ifem_shape ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const unsigned int  i,
const Point p 
)
staticprivate

Definition at line 710 of file fe_interface_inf_fe.C.

References libMesh::INFINITE_MAP, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::LAGRANGE, libMesh::LEGENDRE, libMesh::FEType::radial_family, and libMesh::InfFE< Dim, T_radial, T_map >::shape().

715 {
716  switch (dim)
717  {
718  // 1D
719  case 1:
720  {
721  switch (fe_t.radial_family)
722  {
723  /*
724  * For no derivatives (and local coordinates, as
725  * given in \p p) the infinite element shapes
726  * are independent of mapping type
727  */
728  case INFINITE_MAP:
729  return InfFE<1,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p);
730 
731  case JACOBI_20_00:
732  return InfFE<1,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p);
733 
734  case JACOBI_30_00:
735  return InfFE<1,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p);
736 
737  case LEGENDRE:
738  return InfFE<1,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p);
739 
740  case LAGRANGE:
741  return InfFE<1,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p);
742 
743  default:
744  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
745  }
746  }
747 
748 
749  // 2D
750  case 2:
751  {
752  switch (fe_t.radial_family)
753  {
754  case INFINITE_MAP:
755  return InfFE<2,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p);
756 
757  case JACOBI_20_00:
758  return InfFE<2,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p);
759 
760  case JACOBI_30_00:
761  return InfFE<2,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p);
762 
763  case LEGENDRE:
764  return InfFE<2,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p);
765 
766  case LAGRANGE:
767  return InfFE<2,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p);
768 
769  default:
770  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
771  }
772  }
773 
774 
775  // 3D
776  case 3:
777  {
778  switch (fe_t.radial_family)
779  {
780  case INFINITE_MAP:
781  return InfFE<3,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p);
782 
783  case JACOBI_20_00:
784  return InfFE<3,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p);
785 
786  case JACOBI_30_00:
787  return InfFE<3,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p);
788 
789  case LEGENDRE:
790  return InfFE<3,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p);
791 
792  case LAGRANGE:
793  return InfFE<3,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p);
794 
795  default:
796  libmesh_error_msg("Invalid radial family = " << fe_t.radial_family);
797  }
798  }
799 
800  default:
801  libmesh_error_msg("Invalid dim = " << dim);
802  }
803 }
static Real shape(const FEType &fet, const ElemType t, const unsigned int i, const Point &p)

◆ inverse_map() [1/2]

Point libMesh::FEInterface::inverse_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const Point p,
const Real  tolerance = TOLERANCE,
const bool  secure = true 
)
static
Returns
The location (on the reference element) of the point p located in physical space. This function requires inverting the (probably nonlinear) transformation map, so it is not trivial. The optional parameter tolerance defines how close is "good enough." The map inversion iteration computes the sequence $ \{ p_n \} $, and the iteration is terminated when $ \|p - p_n\| < \mbox{\texttt{tolerance}} $

Definition at line 590 of file fe_interface.C.

References ifem_inverse_map(), is_InfFE_elem(), and libMesh::Elem::type().

Referenced by libMesh::HPCoarsenTest::add_projection(), libMesh::FEMContext::build_new_fe(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_periodic_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_proj_constraints(), libMesh::InfQuad4::contains_point(), libMesh::InfPrism::contains_point(), libMesh::InfHex::contains_point(), libMesh::MeshFunction::discontinuous_gradient(), libMesh::MeshFunction::discontinuous_value(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), inverse_map(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), libMesh::MeshFunction::operator()(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::Elem::point_test(), libMesh::System::point_value(), libMesh::JumpErrorEstimator::reinit_sides(), and libMesh::HPCoarsenTest::select_refinement().

596 {
597 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
598 
599  if (is_InfFE_elem(elem->type()))
600  return ifem_inverse_map(dim, fe_t, elem, p,tolerance, secure);
601 
602 #endif
603 
604  fe_with_vec_switch(inverse_map(elem, p, tolerance, secure));
605 }
static Point ifem_inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
static Point inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_interface.C:590
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ inverse_map() [2/2]

void libMesh::FEInterface::inverse_map ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const std::vector< Point > &  physical_points,
std::vector< Point > &  reference_points,
const Real  tolerance = TOLERANCE,
const bool  secure = true 
)
static
Returns
The location (on the reference element) of the points physical_points located in physical space. This function requires inverting the (probably nonlinear) transformation map, so it is not trivial. The location of each point on the reference element is returned in the vector reference_points. The optional parameter tolerance defines how close is "good enough." The map inversion iteration computes the sequence $ \{ p_n \} $, and the iteration is terminated when $ \|p - p_n\| < \mbox{\texttt{tolerance}} $

Definition at line 610 of file fe_interface.C.

References libMesh::err, ifem_inverse_map(), inverse_map(), is_InfFE_elem(), libMesh::TypeVector< T >::size(), and libMesh::Elem::type().

617 {
618  const std::size_t n_pts = physical_points.size();
619 
620  // Resize the vector
621  reference_points.resize(n_pts);
622 
623  if (n_pts == 0)
624  {
625  libMesh::err << "WARNING: empty vector physical_points!"
626  << std::endl;
627  libmesh_here();
628  return;
629  }
630 
631 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
632 
633  if (is_InfFE_elem(elem->type()))
634  {
635  ifem_inverse_map(dim, fe_t, elem, physical_points, reference_points, tolerance, secure);
636  return;
637  // libmesh_not_implemented();
638  }
639 
640 #endif
641 
642  void_fe_with_vec_switch(inverse_map(elem, physical_points, reference_points, tolerance, secure));
643 }
static Point ifem_inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
static Point inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_interface.C:590
OStreamProxy err(std::cerr)
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ is_InfFE_elem()

bool libMesh::FEInterface::is_InfFE_elem ( const ElemType  et)
staticprivate
Returns
true if et is an element to be processed by class InfFE, false otherwise or when !LIBMESH_ENABLE_INFINITE_ELEMENTS.

Definition at line 45 of file fe_interface.C.

Referenced by compute_data(), inverse_map(), map(), n_dofs(), n_dofs_at_node(), n_dofs_per_elem(), n_shape_functions(), nodal_soln(), and shape().

46 {
47  return false;
48 }

◆ map()

Point libMesh::FEInterface::map ( unsigned int  dim,
const FEType fe_t,
const Elem elem,
const Point p 
)
static
Returns
The point in physical space corresponding to the reference point p which is passed in.

Definition at line 574 of file fe_interface.C.

References ifem_map(), is_InfFE_elem(), and libMesh::Elem::type().

Referenced by libMesh::Elem::point_test().

578 {
579 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
580  if (is_InfFE_elem(elem->type()))
581  return ifem_map(dim, fe_t, elem, p);
582 #endif
583  fe_with_vec_switch(map(elem, p));
584 }
static Point map(unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p)
Definition: fe_interface.C:574
static Point ifem_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p)
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ max_order()

unsigned int libMesh::FEInterface::max_order ( const FEType fe_t,
const ElemType el_t 
)
static
Returns
The maximum polynomial degree that the given finite element family can support on the given geometric element.

Definition at line 1039 of file fe_interface.C.

References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::FEType::family, libMesh::HERMITE, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::HIERARCHIC, libMesh::L2_HIERARCHIC, libMesh::L2_LAGRANGE, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::NEDELEC_ONE, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::QUADSHELL4, libMesh::QUADSHELL8, libMesh::SUBDIVISION, libMesh::SZABAB, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI3SUBDIVISION, libMesh::TRI6, libMesh::TRISHELL3, and libMesh::XYZ.

Referenced by libMesh::DofMap::reinit().

1041 {
1042  // Yeah, I know, infinity is much larger than 11, but our
1043  // solvers don't seem to like high degree polynomials, and our
1044  // quadrature rules and number_lookups tables
1045  // need to go up higher.
1046  const unsigned int unlimited = 11;
1047 
1048  // If we used 0 as a default, then elements missing from this
1049  // table (e.g. infinite elements) would be considered broken.
1050  const unsigned int unknown = unlimited;
1051 
1052  switch (fe_t.family)
1053  {
1054  case LAGRANGE:
1055  case L2_LAGRANGE: // TODO: L2_LAGRANGE can have higher "max_order" than LAGRANGE
1056  case LAGRANGE_VEC:
1057  switch (el_t)
1058  {
1059  case EDGE2:
1060  case EDGE3:
1061  case EDGE4:
1062  return 3;
1063  case TRI3:
1064  case TRISHELL3:
1065  return 1;
1066  case TRI6:
1067  return 2;
1068  case QUAD4:
1069  case QUADSHELL4:
1070  return 1;
1071  case QUAD8:
1072  case QUADSHELL8:
1073  case QUAD9:
1074  return 2;
1075  case TET4:
1076  return 1;
1077  case TET10:
1078  return 2;
1079  case HEX8:
1080  return 1;
1081  case HEX20:
1082  case HEX27:
1083  return 2;
1084  case PRISM6:
1085  return 1;
1086  case PRISM15:
1087  case PRISM18:
1088  return 2;
1089  case PYRAMID5:
1090  return 1;
1091  case PYRAMID13:
1092  case PYRAMID14:
1093  return 2;
1094  default:
1095  return unknown;
1096  }
1097  break;
1098  case MONOMIAL:
1099  switch (el_t)
1100  {
1101  case EDGE2:
1102  case EDGE3:
1103  case EDGE4:
1104  case TRI3:
1105  case TRISHELL3:
1106  case TRI6:
1107  case QUAD4:
1108  case QUADSHELL4:
1109  case QUAD8:
1110  case QUADSHELL8:
1111  case QUAD9:
1112  case TET4:
1113  case TET10:
1114  case HEX8:
1115  case HEX20:
1116  case HEX27:
1117  case PRISM6:
1118  case PRISM15:
1119  case PRISM18:
1120  case PYRAMID5:
1121  case PYRAMID13:
1122  case PYRAMID14:
1123  return unlimited;
1124  default:
1125  return unknown;
1126  }
1127  break;
1128 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
1129  case BERNSTEIN:
1130  switch (el_t)
1131  {
1132  case EDGE2:
1133  case EDGE3:
1134  case EDGE4:
1135  return unlimited;
1136  case TRI3:
1137  case TRISHELL3:
1138  return 1;
1139  case TRI6:
1140  return 6;
1141  case QUAD4:
1142  case QUADSHELL4:
1143  return 1;
1144  case QUAD8:
1145  case QUADSHELL8:
1146  case QUAD9:
1147  return unlimited;
1148  case TET4:
1149  return 1;
1150  case TET10:
1151  return 2;
1152  case HEX8:
1153  return 1;
1154  case HEX20:
1155  return 2;
1156  case HEX27:
1157  return 4;
1158  case PRISM6:
1159  case PRISM15:
1160  case PRISM18:
1161  case PYRAMID5:
1162  case PYRAMID13:
1163  case PYRAMID14:
1164  return 0;
1165  default:
1166  return unknown;
1167  }
1168  break;
1169  case SZABAB:
1170  switch (el_t)
1171  {
1172  case EDGE2:
1173  return 1;
1174  case EDGE3:
1175  case EDGE4:
1176  return 7;
1177  case TRI3:
1178  case TRISHELL3:
1179  return 1;
1180  case TRI6:
1181  return 7;
1182  case QUAD4:
1183  case QUADSHELL4:
1184  return 1;
1185  case QUAD8:
1186  case QUADSHELL8:
1187  case QUAD9:
1188  return 7;
1189  case TET4:
1190  case TET10:
1191  case HEX8:
1192  case HEX20:
1193  case HEX27:
1194  case PRISM6:
1195  case PRISM15:
1196  case PRISM18:
1197  case PYRAMID5:
1198  case PYRAMID13:
1199  case PYRAMID14:
1200  return 0;
1201  default:
1202  return unknown;
1203  }
1204  break;
1205 #endif
1206  case XYZ:
1207  switch (el_t)
1208  {
1209  case EDGE2:
1210  case EDGE3:
1211  case EDGE4:
1212  case TRI3:
1213  case TRISHELL3:
1214  case TRI6:
1215  case QUAD4:
1216  case QUADSHELL4:
1217  case QUAD8:
1218  case QUADSHELL8:
1219  case QUAD9:
1220  case TET4:
1221  case TET10:
1222  case HEX8:
1223  case HEX20:
1224  case HEX27:
1225  case PRISM6:
1226  case PRISM15:
1227  case PRISM18:
1228  case PYRAMID5:
1229  case PYRAMID13:
1230  case PYRAMID14:
1231  return unlimited;
1232  default:
1233  return unknown;
1234  }
1235  break;
1236  case CLOUGH:
1237  switch (el_t)
1238  {
1239  case EDGE2:
1240  case EDGE3:
1241  return 3;
1242  case EDGE4:
1243  case TRI3:
1244  case TRISHELL3:
1245  return 0;
1246  case TRI6:
1247  return 3;
1248  case QUAD4:
1249  case QUADSHELL4:
1250  case QUAD8:
1251  case QUADSHELL8:
1252  case QUAD9:
1253  case TET4:
1254  case TET10:
1255  case HEX8:
1256  case HEX20:
1257  case HEX27:
1258  case PRISM6:
1259  case PRISM15:
1260  case PRISM18:
1261  case PYRAMID5:
1262  case PYRAMID13:
1263  case PYRAMID14:
1264  return 0;
1265  default:
1266  return unknown;
1267  }
1268  break;
1269  case HERMITE:
1270  switch (el_t)
1271  {
1272  case EDGE2:
1273  case EDGE3:
1274  return unlimited;
1275  case EDGE4:
1276  case TRI3:
1277  case TRISHELL3:
1278  case TRI6:
1279  return 0;
1280  case QUAD4:
1281  case QUADSHELL4:
1282  return 3;
1283  case QUAD8:
1284  case QUADSHELL8:
1285  case QUAD9:
1286  return unlimited;
1287  case TET4:
1288  case TET10:
1289  return 0;
1290  case HEX8:
1291  return 3;
1292  case HEX20:
1293  case HEX27:
1294  return unlimited;
1295  case PRISM6:
1296  case PRISM15:
1297  case PRISM18:
1298  case PYRAMID5:
1299  case PYRAMID13:
1300  case PYRAMID14:
1301  return 0;
1302  default:
1303  return unknown;
1304  }
1305  break;
1306  case HIERARCHIC:
1307  switch (el_t)
1308  {
1309  case EDGE2:
1310  case EDGE3:
1311  case EDGE4:
1312  return unlimited;
1313  case TRI3:
1314  case TRISHELL3:
1315  return 1;
1316  case TRI6:
1317  return unlimited;
1318  case QUAD4:
1319  case QUADSHELL4:
1320  return 1;
1321  case QUAD8:
1322  case QUADSHELL8:
1323  case QUAD9:
1324  return unlimited;
1325  case TET4:
1326  case TET10:
1327  return 0;
1328  case HEX8:
1329  case HEX20:
1330  return 1;
1331  case HEX27:
1332  return unlimited;
1333  case PRISM6:
1334  case PRISM15:
1335  case PRISM18:
1336  case PYRAMID5:
1337  case PYRAMID13:
1338  case PYRAMID14:
1339  return 0;
1340  default:
1341  return unknown;
1342  }
1343  break;
1344  case L2_HIERARCHIC:
1345  switch (el_t)
1346  {
1347  case EDGE2:
1348  case EDGE3:
1349  case EDGE4:
1350  return unlimited;
1351  case TRI3:
1352  case TRISHELL3:
1353  return 1;
1354  case TRI6:
1355  return unlimited;
1356  case QUAD4:
1357  case QUADSHELL4:
1358  return 1;
1359  case QUAD8:
1360  case QUADSHELL8:
1361  case QUAD9:
1362  return unlimited;
1363  case TET4:
1364  case TET10:
1365  return 0;
1366  case HEX8:
1367  case HEX20:
1368  return 1;
1369  case HEX27:
1370  return unlimited;
1371  case PRISM6:
1372  case PRISM15:
1373  case PRISM18:
1374  case PYRAMID5:
1375  case PYRAMID13:
1376  case PYRAMID14:
1377  return 0;
1378  default:
1379  return unknown;
1380  }
1381  break;
1382  case SUBDIVISION:
1383  switch (el_t)
1384  {
1385  case TRI3SUBDIVISION:
1386  return unlimited;
1387  default:
1388  return unknown;
1389  }
1390  break;
1391  case NEDELEC_ONE:
1392  switch (el_t)
1393  {
1394  case TRI6:
1395  case QUAD8:
1396  case QUAD9:
1397  case HEX20:
1398  case HEX27:
1399  return 1;
1400  default:
1401  return 0;
1402  }
1403  break;
1404  default:
1405  return 0;
1406  break;
1407  }
1408 }

◆ n_dofs()

unsigned int libMesh::FEInterface::n_dofs ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
static
Returns
The number of shape functions associated with this finite element. Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the total order of the element.

Definition at line 454 of file fe_interface.C.

References ifem_n_dofs(), is_InfFE_elem(), and libMesh::FEType::order.

Referenced by libMesh::DofMap::_dof_indices(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), compute_data(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs(), and libMesh::HPCoarsenTest::select_refinement().

457 {
458 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
459 
460  if (is_InfFE_elem(t))
461  return ifem_n_dofs(dim, fe_t, t);
462 
463 #endif
464 
465  const Order o = fe_t.order;
466 
467  fe_with_vec_switch(n_dofs(t, o));
468 }
static unsigned int n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:454
static unsigned int ifem_n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ n_dofs_at_node()

unsigned int libMesh::FEInterface::n_dofs_at_node ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  n 
)
static
Returns
The number of dofs at node n for a finite element of type fe_t. Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the total order of the element.

Definition at line 473 of file fe_interface.C.

References ifem_n_dofs_at_node(), is_InfFE_elem(), and libMesh::FEType::order.

Referenced by libMesh::DofMap::_dof_indices(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_at_node(), n_dofs_at_node_function(), libMesh::DofMap::old_dof_indices(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), and libMesh::DofMap::reinit().

477 {
478 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
479 
480  if (is_InfFE_elem(t))
481  return ifem_n_dofs_at_node(dim, fe_t, t, n);
482 
483 #endif
484 
485  const Order o = fe_t.order;
486 
487  fe_with_vec_switch(n_dofs_at_node(t, o, n));
488 }
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
Definition: fe_interface.C:473
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45
static unsigned int ifem_n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)

◆ n_dofs_at_node_function()

FEInterface::n_dofs_at_node_ptr libMesh::FEInterface::n_dofs_at_node_function ( const unsigned int  dim,
const FEType fe_t 
)
static
Returns
A function which evaluates n_dofs_at_node for the requested FE type and dimension.

Definition at line 493 of file fe_interface.C.

References n_dofs_at_node().

Referenced by libMesh::DofMap::_dof_indices(), and libMesh::DofMap::old_dof_indices().

495 {
496  fe_with_vec_switch(n_dofs_at_node);
497 }
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
Definition: fe_interface.C:473

◆ n_dofs_per_elem()

unsigned int libMesh::FEInterface::n_dofs_per_elem ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
static
Returns
The number of dofs interior to the element, not associated with any interior nodes. Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the total order of the element.

Definition at line 504 of file fe_interface.C.

References ifem_n_dofs_per_elem(), is_InfFE_elem(), and libMesh::FEType::order.

Referenced by libMesh::DofMap::_dof_indices(), libMesh::InfFE< Dim, T_radial, T_map >::compute_shape_indices(), libMesh::InfFE< Dim, T_radial, T_map >::n_dofs_per_elem(), libMesh::DofMap::old_dof_indices(), and libMesh::DofMap::reinit().

507 {
508 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
509 
510  if (is_InfFE_elem(t))
511  return ifem_n_dofs_per_elem(dim, fe_t, t);
512 
513 #endif
514 
515  const Order o = fe_t.order;
516 
517  fe_with_vec_switch(n_dofs_per_elem(t, o));
518 }
static unsigned int n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:504
static unsigned int ifem_n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ n_shape_functions()

unsigned int libMesh::FEInterface::n_shape_functions ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t 
)
static
Returns
The number of shape functions associated with this finite element of type fe_t. Automatically decides which finite element class to use.

On a p-refined element, fe_t.order should be the total order of the element.

Definition at line 428 of file fe_interface.C.

References ifem_n_shape_functions(), is_InfFE_elem(), and libMesh::FEType::order.

431 {
432 
433 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
434  /*
435  * Since the FEType, stored in DofMap/(some System child), has to
436  * be the _same_ for InfFE and FE, we have to catch calls
437  * to infinite elements through the element type.
438  */
439 
440  if (is_InfFE_elem(t))
441  return ifem_n_shape_functions(dim, fe_t, t);
442 
443 #endif
444 
445  const Order o = fe_t.order;
446 
447  fe_with_vec_switch(n_shape_functions(t, o));
448 }
static unsigned int ifem_n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:428
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ n_vec_dim()

unsigned int libMesh::FEInterface::n_vec_dim ( const MeshBase mesh,
const FEType fe_type 
)
static
Returns
The number of components of a vector-valued element. Scalar-valued elements return 1.

Definition at line 1454 of file fe_interface.C.

References libMesh::FEType::family, libMesh::LAGRANGE_VEC, mesh, and libMesh::NEDELEC_ONE.

Referenced by libMesh::ExactSolution::_compute_error(), libMesh::EquationSystems::build_parallel_solution_vector(), and libMesh::EquationSystems::build_variable_names().

1456 {
1457  switch (fe_type.family)
1458  {
1459  //FIXME: We currently assume that the number of vector components is tied
1460  // to the mesh dimension. This will break for mixed-dimension meshes.
1461  case LAGRANGE_VEC:
1462  case NEDELEC_ONE:
1463  return mesh.mesh_dimension();
1464  default:
1465  return 1;
1466  }
1467 }
MeshBase & mesh

◆ nodal_soln()

void libMesh::FEInterface::nodal_soln ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const std::vector< Number > &  elem_soln,
std::vector< Number > &  nodal_soln 
)
static

Build the nodal soln from the element soln. This is the solution that will be plotted. Automatically passes the request to the appropriate finite element class member. To indicate that results from this specific implementation of nodal_soln should not be used, the vector nodal_soln is returned empty.

On a p-refined element, fe_t.order should be the base order of the element.

Definition at line 550 of file fe_interface.C.

References ifem_nodal_soln(), is_InfFE_elem(), libMesh::FEType::order, and libMesh::Elem::type().

Referenced by libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), ifem_nodal_soln(), libMesh::EnsightIO::write_scalar_ascii(), and libMesh::EnsightIO::write_vector_ascii().

555 {
556 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
557 
558  if (is_InfFE_elem(elem->type()))
559  {
560  ifem_nodal_soln(dim, fe_t, elem, elem_soln, nodal_soln);
561  return;
562  }
563 
564 #endif
565 
566  const Order order = fe_t.order;
567 
568  void_fe_with_vec_switch(nodal_soln(elem, order, elem_soln, nodal_soln));
569 }
static void ifem_nodal_soln(const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45
static void nodal_soln(const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
Definition: fe_interface.C:550

◆ on_reference_element()

bool libMesh::FEInterface::on_reference_element ( const Point p,
const ElemType  t,
const Real  eps = TOLERANCE 
)
static
Returns
true if the point p is located on the reference element for element type t, false otherwise.

Since we are doing floating point comparisons, the parameter eps can be specified to indicate a tolerance. For example, $ \xi \le 1 $ becomes $ \xi \le 1 + \epsilon $.

Definition at line 647 of file fe_interface.C.

References libMesh::FEAbstract::on_reference_element().

Referenced by libMesh::InfQuad4::contains_point(), libMesh::InfPrism::contains_point(), libMesh::InfHex::contains_point(), and libMesh::Elem::point_test().

650 {
651  return FEBase::on_reference_element(p,t,eps);
652 }
static bool on_reference_element(const Point &p, const ElemType t, const Real eps=TOLERANCE)
Definition: fe_abstract.C:581

◆ shape() [1/8]

Real libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  i,
const Point p 
)
static
Returns
The value of the $ i^{th} $ shape function at point p. This method allows you to specify the dimension, element type, and order directly. Automatically passes the request to the appropriate finite element class member.
Note
On a p-refined element, fe_t.order should be the total order of the element.

Definition at line 657 of file fe_interface.C.

References ifem_shape(), is_InfFE_elem(), and libMesh::FEType::order.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::compute_data(), compute_data(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), shape(), and libMesh::InfFE< Dim, T_radial, T_map >::shape().

662 {
663 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
664 
665  if (is_InfFE_elem(t))
666  return ifem_shape(dim, fe_t, t, i, p);
667 
668 #endif
669 
670  const Order o = fe_t.order;
671 
672  fe_switch(shape(t,o,i,p));
673 }
static Real ifem_shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
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:657
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ shape() [2/8]

Real libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const unsigned int  i,
const Point p 
)
static
Returns
The value of the $ i^{th} $ shape function at point p. This method allows you to specify the dimension, element type, and order directly. Automatically passes the request to the appropriate finite element class member.
Note
On a p-refined element, fe_t.order should be the base order of the element.

Definition at line 675 of file fe_interface.C.

References ifem_shape(), is_InfFE_elem(), libMesh::FEType::order, shape(), and libMesh::Elem::type().

680 {
681 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
682 
683  if (elem && is_InfFE_elem(elem->type()))
684  return ifem_shape(dim, fe_t, elem, i, p);
685 
686 #endif
687 
688  const Order o = fe_t.order;
689 
690  fe_switch(shape(elem,o,i,p));
691 }
static Real ifem_shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
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:657
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ shape() [3/8]

template<typename OutputType >
static void libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  i,
const Point p,
OutputType &  phi 
)
static
Returns
The value of the $ i^{th} $ shape function at point p. This method allows you to specify the dimension, element type, and order directly. Automatically passes the request to the appropriate scalar finite element class member.
Note
On a p-refined element, fe_t.order should be the total order of the element.

◆ shape() [4/8]

template<typename OutputType >
static void libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const unsigned int  i,
const Point p,
OutputType &  phi 
)
static
Returns
The value of the $ i^{th} $ shape function at point p. This method allows you to specify the dimension, element type, and order directly. Automatically passes the request to the appropriate scalar finite element class member.
Note
On a p-refined element, fe_t.order should be the total order of the element.

◆ shape() [5/8]

template<>
void libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  i,
const Point p,
Real phi 
)

Definition at line 694 of file fe_interface.C.

700 {
701 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
702 
703  if (is_InfFE_elem(t))
704  {
705  phi = ifem_shape(dim, fe_t, t, i, p);
706  return;
707  }
708 
709 #endif
710 
711  const Order o = fe_t.order;
712 
713  switch(dim)
714  {
715  case 0:
716  fe_scalar_vec_error_switch(0, shape(t,o,i,p), phi = , ; break;);
717  break;
718  case 1:
719  fe_scalar_vec_error_switch(1, shape(t,o,i,p), phi = , ; break;);
720  break;
721  case 2:
722  fe_scalar_vec_error_switch(2, shape(t,o,i,p), phi = , ; break;);
723  break;
724  case 3:
725  fe_scalar_vec_error_switch(3, shape(t,o,i,p), phi = , ; break;);
726  break;
727  default:
728  libmesh_error_msg("Invalid dimension = " << dim);
729  }
730 
731  return;
732 }
static Real ifem_shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
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:657
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ shape() [6/8]

template<>
void libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const unsigned int  i,
const Point p,
Real phi 
)

Definition at line 735 of file fe_interface.C.

741 {
742 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
743 
744  if (elem && is_InfFE_elem(elem->type()))
745  {
746  phi = ifem_shape(dim, fe_t, elem, i, p);
747  return;
748  }
749 #endif
750 
751  const Order o = fe_t.order;
752 
753  switch(dim)
754  {
755  case 0:
756  fe_scalar_vec_error_switch(0, shape(elem,o,i,p), phi = , ; break;);
757  break;
758  case 1:
759  fe_scalar_vec_error_switch(1, shape(elem,o,i,p), phi = , ; break;);
760  break;
761  case 2:
762  fe_scalar_vec_error_switch(2, shape(elem,o,i,p), phi = , ; break;);
763  break;
764  case 3:
765  fe_scalar_vec_error_switch(3, shape(elem,o,i,p), phi = , ; break;);
766  break;
767  default:
768  libmesh_error_msg("Invalid dimension = " << dim);
769  }
770 
771  return;
772 }
static Real ifem_shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
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:657
static bool is_InfFE_elem(const ElemType et)
Definition: fe_interface.C:45

◆ shape() [7/8]

template<>
void libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const ElemType  t,
const unsigned int  i,
const Point p,
RealGradient phi 
)

Definition at line 775 of file fe_interface.C.

781 {
782  // This even does not handle infinite elements at all!?
783  const Order o = fe_t.order;
784 
785  switch(dim)
786  {
787  case 0:
788  fe_vector_scalar_error_switch(0, shape(t,o,i,p), phi = , ; break;);
789  break;
790  case 1:
791  fe_vector_scalar_error_switch(1, shape(t,o,i,p), phi = , ; break;);
792  break;
793  case 2:
794  fe_vector_scalar_error_switch(2, shape(t,o,i,p), phi = , ; break;);
795  break;
796  case 3:
797  fe_vector_scalar_error_switch(3, shape(t,o,i,p), phi = , ; break;);
798  break;
799  default:
800  libmesh_error_msg("Invalid dimension = " << dim);
801  }
802 
803  return;
804 }
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:657

◆ shape() [8/8]

template<>
void libMesh::FEInterface::shape ( const unsigned int  dim,
const FEType fe_t,
const Elem elem,
const unsigned int  i,
const Point p,
RealGradient phi 
)

Definition at line 807 of file fe_interface.C.

813 {
814  const Order o = fe_t.order;
815 
816  switch(dim)
817  {
818  case 0:
819  fe_vector_scalar_error_switch(0, shape(elem,o,i,p), phi = , ; break;);
820  break;
821  case 1:
822  fe_vector_scalar_error_switch(1, shape(elem,o,i,p), phi = , ; break;);
823  break;
824  case 2:
825  fe_vector_scalar_error_switch(2, shape(elem,o,i,p), phi = , ; break;);
826  break;
827  case 3:
828  fe_vector_scalar_error_switch(3, shape(elem,o,i,p), phi = , ; break;);
829  break;
830  default:
831  libmesh_error_msg("Invalid dimension = " << dim);
832  }
833 
834  return;
835 }
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:657

The documentation for this class was generated from the following files: