libMesh::FEInterface Class Reference
#include <fe_interface.h>
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 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_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
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.
Definition at line 64 of file fe_interface.h.
Constructor & Destructor Documentation
| libMesh::FEInterface::FEInterface | ( | ) | [private] |
Empty constructor. Do not create an object of this type.
Definition at line 32 of file fe_interface.C.
References libMesh::err.
00033 { 00034 libMesh::err << "ERROR: Do not define an object of this type." 00035 << std::endl; 00036 libmesh_error(); 00037 }
| virtual libMesh::FEInterface::~FEInterface | ( | ) | [inline, virtual] |
Member Function Documentation
| 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 826 of file fe_interface.C.
References libMeshEnums::CLOUGH, libMesh::Elem::dim(), libMesh::FEType::family, libMeshEnums::HERMITE, libMeshEnums::HIERARCHIC, libMeshEnums::L2_HIERARCHIC, libMeshEnums::LAGRANGE, libMeshEnums::LAGRANGE_VEC, and libMesh::DofMap::variable_type().
00830 { 00831 libmesh_assert(elem); 00832 00833 const FEType& fe_t = dof_map.variable_type(variable_number); 00834 00835 switch (elem->dim()) 00836 { 00837 case 0: 00838 case 1: 00839 { 00840 // No constraints in 0D/1D. 00841 return; 00842 } 00843 00844 00845 case 2: 00846 { 00847 switch (fe_t.family) 00848 { 00849 case CLOUGH: 00850 FE<2,CLOUGH>::compute_constraints (constraints, 00851 dof_map, 00852 variable_number, 00853 elem); return; 00854 00855 case HERMITE: 00856 FE<2,HERMITE>::compute_constraints (constraints, 00857 dof_map, 00858 variable_number, 00859 elem); return; 00860 00861 case LAGRANGE: 00862 FE<2,LAGRANGE>::compute_constraints (constraints, 00863 dof_map, 00864 variable_number, 00865 elem); return; 00866 00867 case HIERARCHIC: 00868 FE<2,HIERARCHIC>::compute_constraints (constraints, 00869 dof_map, 00870 variable_number, 00871 elem); return; 00872 00873 case L2_HIERARCHIC: 00874 FE<2,L2_HIERARCHIC>::compute_constraints (constraints, 00875 dof_map, 00876 variable_number, 00877 elem); return; 00878 00879 case LAGRANGE_VEC: 00880 FE<2,LAGRANGE_VEC>::compute_constraints (constraints, 00881 dof_map, 00882 variable_number, 00883 elem); return; 00884 00885 00886 default: 00887 return; 00888 } 00889 } 00890 00891 00892 case 3: 00893 { 00894 switch (fe_t.family) 00895 { 00896 case HERMITE: 00897 FE<3,HERMITE>::compute_constraints (constraints, 00898 dof_map, 00899 variable_number, 00900 elem); return; 00901 00902 case LAGRANGE: 00903 FE<3,LAGRANGE>::compute_constraints (constraints, 00904 dof_map, 00905 variable_number, 00906 elem); return; 00907 00908 case HIERARCHIC: 00909 FE<3,HIERARCHIC>::compute_constraints (constraints, 00910 dof_map, 00911 variable_number, 00912 elem); return; 00913 00914 case L2_HIERARCHIC: 00915 FE<3,L2_HIERARCHIC>::compute_constraints (constraints, 00916 dof_map, 00917 variable_number, 00918 elem); return; 00919 00920 case LAGRANGE_VEC: 00921 FE<3,LAGRANGE_VEC>::compute_constraints (constraints, 00922 dof_map, 00923 variable_number, 00924 elem); return; 00925 default: 00926 return; 00927 } 00928 } 00929 00930 00931 default: 00932 libmesh_error(); 00933 } 00934 }
| 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 returns the values also through data. See this as a generalization of shape(). Currently, with disabled infinite elements, returns a vector of all shape functions of elem evaluated ap p.
On a p-refined element, fe_t.order should be the base order of the element.
Definition at line 790 of file fe_interface.C.
References ifem_compute_data(), libMesh::FEComputeData::init(), is_InfFE_elem(), n_dofs(), libMesh::FEType::order, libMesh::FEComputeData::p, libMesh::Elem::p_level(), shape(), libMesh::FEComputeData::shape, and libMesh::Elem::type().
Referenced by ifem_compute_data(), and libMesh::MeshFunction::operator()().
00794 { 00795 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00796 00797 if ( is_InfFE_elem(elem->type()) ) 00798 { 00799 data.init(); 00800 ifem_compute_data(dim, fe_t, elem, data); 00801 return; 00802 } 00803 00804 #endif 00805 00806 FEType p_refined = fe_t; 00807 p_refined.order = static_cast<Order>(p_refined.order + elem->p_level()); 00808 00809 const unsigned int n_dof = n_dofs (dim, p_refined, elem->type()); 00810 const Point& p = data.p; 00811 data.shape.resize(n_dof); 00812 00813 // set default values for all the output fields 00814 data.init(); 00815 00816 for (unsigned int n=0; n<n_dof; n++) 00817 data.shape[n] = shape(dim, p_refined, elem, n, p); 00818 00819 return; 00820 }
| 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 942 of file fe_interface.C.
00949 { 00950 // No element-specific optimizations currently exist 00951 FEBase::compute_periodic_constraints (constraints, 00952 dof_map, 00953 boundaries, 00954 mesh, 00955 point_locator, 00956 variable_number, 00957 elem); 00958 }
| 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 487 of file fe_interface.C.
References libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::BoundaryProjectSolution::operator()(), libMesh::ProjectFEMSolution::operator()(), and libMesh::ProjectSolution::operator()().
00492 { 00493 const Order o = fe_t.order; 00494 00495 void_fe_with_vec_switch(dofs_on_edge(elem, o, e, di)); 00496 00497 libmesh_error(); 00498 }
| 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 472 of file fe_interface.C.
References libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), libMesh::FEGenericBase< OutputType >::compute_proj_constraints(), libMesh::BoundaryProjectSolution::operator()(), libMesh::ProjectFEMSolution::operator()(), and libMesh::ProjectSolution::operator()().
00477 { 00478 const Order o = fe_t.order; 00479 00480 void_fe_with_vec_switch(dofs_on_side(elem, o, s, di)); 00481 00482 libmesh_error(); 00483 }
| 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 1279 of file fe_interface.C.
References libMeshEnums::BERNSTEIN, libMeshEnums::CLOUGH, libMesh::FEType::family, libMeshEnums::HERMITE, libMeshEnums::HIERARCHIC, libMeshEnums::L2_HIERARCHIC, libMeshEnums::L2_LAGRANGE, libMeshEnums::LAGRANGE, libMeshEnums::LAGRANGE_VEC, libMeshEnums::MONOMIAL, libMeshEnums::NEDELEC_ONE, libMeshEnums::SZABAB, and libMeshEnums::XYZ.
Referenced by libMesh::DofMap::dof_indices(), libMesh::DofMap::old_dof_indices(), and libMesh::DofMap::reinit().
01280 { 01281 switch (fe_t.family) 01282 { 01283 case LAGRANGE: 01284 case L2_LAGRANGE: 01285 case MONOMIAL: 01286 case L2_HIERARCHIC: 01287 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 01288 case BERNSTEIN: 01289 case SZABAB: 01290 #endif 01291 case XYZ: 01292 case LAGRANGE_VEC: 01293 case NEDELEC_ONE: 01294 return false; 01295 case CLOUGH: 01296 case HERMITE: 01297 case HIERARCHIC: 01298 default: 01299 return true; 01300 } 01301 }
| 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 1308 of file fe_interface.C.
References libMeshEnums::LAGRANGE_VEC, libMeshEnums::NEDELEC_ONE, libMeshEnums::TYPE_SCALAR, and libMeshEnums::TYPE_VECTOR.
01309 { 01310 switch (fe_family) 01311 { 01312 case LAGRANGE_VEC: 01313 case NEDELEC_ONE: 01314 return TYPE_VECTOR; 01315 default: 01316 return TYPE_SCALAR; 01317 } 01318 }
| FEFieldType libMesh::FEInterface::field_type | ( | const FEType & | fe_type | ) | [static] |
Returns the number of components of a vector-valued element. Scalar-valued elements return 1.
Definition at line 1303 of file fe_interface.C.
References libMesh::FEType::family.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::EquationSystems::build_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::ExactSolution::compute_error(), libMesh::ExactSolution::error_norm(), libMesh::FEMContext::FEMContext(), and libMesh::FE< Dim, T >::init_shape_functions().
01304 { 01305 return FEInterface::field_type( fe_type.family ); 01306 }
| void libMesh::FEInterface::ifem_compute_data | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const Elem * | elem, | |||
| FEComputeData & | data | |||
| ) | [static, private] |
Definition at line 871 of file fe_interface_inf_fe.C.
References compute_data(), libMeshEnums::INFINITE_MAP, libMeshEnums::JACOBI_20_00, libMeshEnums::JACOBI_30_00, libMeshEnums::LAGRANGE, libMeshEnums::LEGENDRE, and libMesh::FEType::radial_family.
Referenced by compute_data().
00875 { 00876 switch (dim) 00877 { 00878 // 1D 00879 case 1: 00880 { 00881 switch (fe_t.radial_family) 00882 { 00883 /* 00884 * For no derivatives (and local coordinates, as 00885 * given in \p p) the infinite element shapes 00886 * are independent of mapping type 00887 */ 00888 case INFINITE_MAP: 00889 InfFE<1,INFINITE_MAP,CARTESIAN>::compute_data(fe_t, elem, data); 00890 break; 00891 00892 case JACOBI_20_00: 00893 InfFE<1,JACOBI_20_00,CARTESIAN>::compute_data(fe_t, elem, data); 00894 break; 00895 00896 case JACOBI_30_00: 00897 InfFE<1,JACOBI_30_00,CARTESIAN>::compute_data(fe_t, elem, data); 00898 break; 00899 00900 case LEGENDRE: 00901 InfFE<1,LEGENDRE,CARTESIAN>::compute_data(fe_t, elem, data); 00902 break; 00903 00904 case LAGRANGE: 00905 InfFE<1,LAGRANGE,CARTESIAN>::compute_data(fe_t, elem, data); 00906 break; 00907 00908 default: 00909 libmesh_error(); 00910 } 00911 00912 break; 00913 } 00914 00915 00916 // 2D 00917 case 2: 00918 { 00919 switch (fe_t.radial_family) 00920 { 00921 case INFINITE_MAP: 00922 InfFE<2,INFINITE_MAP,CARTESIAN>::compute_data(fe_t, elem, data); 00923 break; 00924 00925 case JACOBI_20_00: 00926 InfFE<2,JACOBI_20_00,CARTESIAN>::compute_data(fe_t, elem, data); 00927 break; 00928 00929 case JACOBI_30_00: 00930 InfFE<2,JACOBI_30_00,CARTESIAN>::compute_data(fe_t, elem, data); 00931 break; 00932 00933 case LEGENDRE: 00934 InfFE<2,LEGENDRE,CARTESIAN>::compute_data(fe_t, elem, data); 00935 break; 00936 00937 case LAGRANGE: 00938 InfFE<2,LAGRANGE,CARTESIAN>::compute_data(fe_t, elem, data); 00939 break; 00940 00941 default: 00942 libmesh_error(); 00943 } 00944 00945 break; 00946 } 00947 00948 00949 // 3D 00950 case 3: 00951 { 00952 switch (fe_t.radial_family) 00953 { 00954 case INFINITE_MAP: 00955 InfFE<3,INFINITE_MAP,CARTESIAN>::compute_data(fe_t, elem, data); 00956 break; 00957 00958 case JACOBI_20_00: 00959 InfFE<3,JACOBI_20_00,CARTESIAN>::compute_data(fe_t, elem, data); 00960 break; 00961 00962 case JACOBI_30_00: 00963 InfFE<3,JACOBI_30_00,CARTESIAN>::compute_data(fe_t, elem, data); 00964 break; 00965 00966 case LEGENDRE: 00967 InfFE<3,LEGENDRE,CARTESIAN>::compute_data(fe_t, elem, data); 00968 break; 00969 00970 case LAGRANGE: 00971 InfFE<3,LAGRANGE,CARTESIAN>::compute_data(fe_t, elem, data); 00972 break; 00973 00974 default: 00975 libmesh_error(); 00976 } 00977 00978 break; 00979 } 00980 00981 00982 default: 00983 libmesh_error(); 00984 break; 00985 } 00986 00987 return; 00988 }
| 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 | |||
| ) | [static, private] |
Definition at line 583 of file fe_interface_inf_fe.C.
References libMeshEnums::CARTESIAN, libMesh::FEType::inf_map, and inverse_map().
00590 { 00591 switch (dim) 00592 { 00593 // 1D 00594 case 1: 00595 { 00596 switch (fe_t.inf_map) 00597 { 00598 case CARTESIAN: 00599 InfFE<1,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure); 00600 return; 00601 00602 default: 00603 libmesh_error(); 00604 } 00605 } 00606 00607 00608 // 2D 00609 case 2: 00610 { 00611 switch (fe_t.inf_map) 00612 { 00613 case CARTESIAN: 00614 InfFE<2,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure); 00615 return; 00616 00617 default: 00618 libmesh_error(); 00619 } 00620 00621 } 00622 00623 00624 // 3D 00625 case 3: 00626 { 00627 switch (fe_t.inf_map) 00628 { 00629 case CARTESIAN: 00630 InfFE<3,JACOBI_20_00,CARTESIAN>::inverse_map(elem, physical_points, reference_points, tolerance, secure); 00631 return; 00632 00633 default: 00634 libmesh_error(); 00635 } 00636 00637 } 00638 00639 00640 default: 00641 libmesh_error(); 00642 } 00643 00644 libmesh_error(); 00645 return; 00646 }
| 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 | |||
| ) | [static, private] |
Definition at line 470 of file fe_interface_inf_fe.C.
References libMeshEnums::CARTESIAN, libMeshEnums::ELLIPSOIDAL, libMesh::err, libMesh::FEType::inf_map, inverse_map(), and libMeshEnums::SPHERICAL.
Referenced by inverse_map().
00476 { 00477 switch (dim) 00478 { 00479 // 1D 00480 case 1: 00481 { 00482 switch (fe_t.inf_map) 00483 { 00484 case CARTESIAN: 00485 return InfFE<1,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure); 00486 00487 case SPHERICAL: 00488 case ELLIPSOIDAL: 00489 { 00490 libMesh::err << "ERROR: Spherical and Ellipsoidal IFEMs not (yet) " << std::endl 00491 << "implemented." << std::endl; 00492 libmesh_error(); 00493 } 00494 00495 /* 00496 case SPHERICAL: 00497 return InfFE<1,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance); 00498 00499 case ELLIPSOIDAL: 00500 return InfFE<1,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance); 00501 */ 00502 00503 default: 00504 libmesh_error(); 00505 } 00506 } 00507 00508 00509 // 2D 00510 case 2: 00511 { 00512 switch (fe_t.inf_map) 00513 { 00514 case CARTESIAN: 00515 return InfFE<2,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure); 00516 00517 case SPHERICAL: 00518 case ELLIPSOIDAL: 00519 { 00520 libMesh::err << "ERROR: Spherical and Ellipsoidal IFEMs not (yet) " << std::endl 00521 << "implemented." << std::endl; 00522 libmesh_error(); 00523 } 00524 00525 /* 00526 case SPHERICAL: 00527 return InfFE<2,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance); 00528 00529 case ELLIPSOIDAL: 00530 return InfFE<2,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance); 00531 */ 00532 00533 default: 00534 libmesh_error(); 00535 } 00536 00537 } 00538 00539 00540 // 3D 00541 case 3: 00542 { 00543 switch (fe_t.inf_map) 00544 { 00545 case CARTESIAN: 00546 return InfFE<3,JACOBI_20_00,CARTESIAN>::inverse_map(elem, p, tolerance, secure); 00547 00548 case SPHERICAL: 00549 case ELLIPSOIDAL: 00550 { 00551 libMesh::err << "ERROR: Spherical and Ellipsoidal IFEMs not (yet) " << std::endl 00552 << "implemented." << std::endl; 00553 libmesh_error(); 00554 } 00555 00556 /* 00557 case SPHERICAL: 00558 return InfFE<3,JACOBI_20_00,SPHERICAL>::inverse_map(elem, p, tolerance); 00559 00560 case ELLIPSOIDAL: 00561 return InfFE<3,JACOBI_20_00,ELLIPSOIDAL>::inverse_map(elem, p, tolerance); 00562 */ 00563 00564 default: 00565 libmesh_error(); 00566 } 00567 00568 } 00569 00570 00571 default: 00572 libmesh_error(); 00573 } 00574 00575 00576 libmesh_error(); 00577 Point pt; 00578 return pt; 00579 }
| unsigned int libMesh::FEInterface::ifem_n_dofs | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const ElemType | t | |||
| ) | [static, private] |
Definition at line 75 of file fe_interface_inf_fe.C.
References n_dofs().
Referenced by n_dofs().
00078 { 00079 switch (dim) 00080 { 00081 // 1D 00082 case 1: 00083 /* 00084 * Since InfFE<Dim,T_radial,T_map>::n_dofs(...) 00085 * is actually independent of T_radial and T_map, we can use 00086 * just any T_radial and T_map 00087 */ 00088 return InfFE<1,JACOBI_20_00,CARTESIAN>::n_dofs(fe_t, t); 00089 00090 // 2D 00091 case 2: 00092 return InfFE<2,JACOBI_20_00,CARTESIAN>::n_dofs(fe_t, t); 00093 00094 // 3D 00095 case 3: 00096 return InfFE<3,JACOBI_20_00,CARTESIAN>::n_dofs(fe_t, t); 00097 00098 default: 00099 libmesh_error(); 00100 } 00101 00102 00103 libmesh_error(); 00104 return 0; 00105 }
| unsigned int libMesh::FEInterface::ifem_n_dofs_at_node | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const ElemType | t, | |||
| const unsigned int | n | |||
| ) | [static, private] |
Definition at line 110 of file fe_interface_inf_fe.C.
References n_dofs_at_node().
Referenced by n_dofs_at_node().
00114 { 00115 switch (dim) 00116 { 00117 // 1D 00118 case 1: 00119 /* 00120 * Since InfFE<Dim,T_radial,T_map>::n_dofs_at_node(...) 00121 * is actually independent of T_radial and T_map, we can use 00122 * just any T_radial and T_map 00123 */ 00124 return InfFE<1,JACOBI_20_00,CARTESIAN>::n_dofs_at_node(fe_t, t, n); 00125 00126 // 2D 00127 case 2: 00128 return InfFE<2,JACOBI_20_00,CARTESIAN>::n_dofs_at_node(fe_t, t, n); 00129 00130 // 3D 00131 case 3: 00132 return InfFE<3,JACOBI_20_00,CARTESIAN>::n_dofs_at_node(fe_t, t, n); 00133 00134 default: 00135 libmesh_error(); 00136 } 00137 00138 00139 libmesh_error(); 00140 return 0; 00141 }
| unsigned int libMesh::FEInterface::ifem_n_dofs_per_elem | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const ElemType | t | |||
| ) | [static, private] |
Definition at line 147 of file fe_interface_inf_fe.C.
References n_dofs_per_elem().
Referenced by n_dofs_per_elem().
00150 { 00151 switch (dim) 00152 { 00153 // 1D 00154 case 1: 00155 /* 00156 * Since InfFE<Dim,T_radial,T_map>::n_dofs(...) 00157 * is actually independent of T_radial and T_map, we can use 00158 * just any T_radial and T_map 00159 */ 00160 return InfFE<1,JACOBI_20_00,CARTESIAN>::n_dofs_per_elem(fe_t, t); 00161 00162 // 2D 00163 case 2: 00164 return InfFE<2,JACOBI_20_00,CARTESIAN>::n_dofs_per_elem(fe_t, t); 00165 00166 // 3D 00167 case 3: 00168 return InfFE<3,JACOBI_20_00,CARTESIAN>::n_dofs_per_elem(fe_t, t); 00169 00170 default: 00171 libmesh_error(); 00172 } 00173 00174 00175 libmesh_error(); 00176 return 0; 00177 }
| unsigned int libMesh::FEInterface::ifem_n_shape_functions | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const ElemType | t | |||
| ) | [static, private] |
Definition at line 39 of file fe_interface_inf_fe.C.
References n_shape_functions().
Referenced by n_shape_functions().
00042 { 00043 switch (dim) 00044 { 00045 // 1D 00046 case 1: 00047 /* 00048 * Since InfFE<Dim,T_radial,T_map>::n_shape_functions(...) 00049 * is actually independent of T_radial and T_map, we can use 00050 * just any T_radial and T_map 00051 */ 00052 return InfFE<1,JACOBI_20_00,CARTESIAN>::n_shape_functions(fe_t, t); 00053 00054 // 2D 00055 case 2: 00056 return InfFE<2,JACOBI_20_00,CARTESIAN>::n_shape_functions(fe_t, t); 00057 00058 // 3D 00059 case 3: 00060 return InfFE<3,JACOBI_20_00,CARTESIAN>::n_shape_functions(fe_t, t); 00061 00062 default: 00063 libmesh_error(); 00064 } 00065 00066 00067 libmesh_error(); 00068 return 0; 00069 }
| 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 | |||
| ) | [static, private] |
Definition at line 182 of file fe_interface_inf_fe.C.
References libMeshEnums::CARTESIAN, libMesh::err, libMesh::FEType::inf_map, libMeshEnums::INFINITE_MAP, libMeshEnums::JACOBI_20_00, libMeshEnums::JACOBI_30_00, libMeshEnums::LAGRANGE, libMeshEnums::LEGENDRE, nodal_soln(), and libMesh::FEType::radial_family.
Referenced by nodal_soln().
00187 { 00188 switch (dim) 00189 { 00190 00191 // 1D 00192 case 1: 00193 { 00194 switch (fe_t.radial_family) 00195 { 00196 case INFINITE_MAP: 00197 { 00198 libMesh::err << "ERROR: INFINTE_MAP is not a valid shape family for radial approximation." << std::endl; 00199 libmesh_error(); 00200 break; 00201 } 00202 00203 case JACOBI_20_00: 00204 { 00205 switch (fe_t.inf_map) 00206 { 00207 case CARTESIAN: 00208 { 00209 InfFE<1,JACOBI_20_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00210 break; 00211 } 00212 default: 00213 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00214 libmesh_error(); 00215 } 00216 break; 00217 } 00218 00219 case JACOBI_30_00: 00220 { 00221 switch (fe_t.inf_map) 00222 { 00223 case CARTESIAN: 00224 { 00225 InfFE<1,JACOBI_30_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00226 break; 00227 } 00228 default: 00229 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00230 libmesh_error(); 00231 } 00232 break; 00233 } 00234 00235 case LEGENDRE: 00236 { 00237 switch (fe_t.inf_map) 00238 { 00239 case CARTESIAN: 00240 { 00241 InfFE<1,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00242 break; 00243 } 00244 default: 00245 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00246 libmesh_error(); 00247 } 00248 break; 00249 } 00250 00251 case LAGRANGE: 00252 { 00253 switch (fe_t.inf_map) 00254 { 00255 case CARTESIAN: 00256 { 00257 InfFE<1,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00258 break; 00259 } 00260 default: 00261 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00262 libmesh_error(); 00263 } 00264 break; 00265 } 00266 00267 00268 00269 default: 00270 libMesh::err << "ERROR: Bad FEType.radial_family= " << fe_t.radial_family << std::endl; 00271 libmesh_error(); 00272 break; 00273 } 00274 00275 break; 00276 } 00277 00278 00279 00280 00281 // 2D 00282 case 2: 00283 { 00284 switch (fe_t.radial_family) 00285 { 00286 case INFINITE_MAP: 00287 { 00288 libMesh::err << "ERROR: INFINTE_MAP is not a valid shape family for radial approximation." << std::endl; 00289 libmesh_error(); 00290 break; 00291 } 00292 00293 case JACOBI_20_00: 00294 { 00295 switch (fe_t.inf_map) 00296 { 00297 case CARTESIAN: 00298 { 00299 InfFE<2,JACOBI_20_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00300 break; 00301 } 00302 default: 00303 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00304 libmesh_error(); 00305 } 00306 break; 00307 } 00308 00309 case JACOBI_30_00: 00310 { 00311 switch (fe_t.inf_map) 00312 { 00313 case CARTESIAN: 00314 { 00315 InfFE<2,JACOBI_30_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00316 break; 00317 } 00318 default: 00319 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00320 libmesh_error(); 00321 } 00322 break; 00323 } 00324 00325 case LEGENDRE: 00326 { 00327 switch (fe_t.inf_map) 00328 { 00329 case CARTESIAN: 00330 { 00331 InfFE<2,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00332 break; 00333 } 00334 default: 00335 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00336 libmesh_error(); 00337 } 00338 break; 00339 } 00340 00341 case LAGRANGE: 00342 { 00343 switch (fe_t.inf_map) 00344 { 00345 case CARTESIAN: 00346 { 00347 InfFE<2,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00348 break; 00349 } 00350 default: 00351 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00352 libmesh_error(); 00353 } 00354 break; 00355 } 00356 00357 00358 00359 default: 00360 libMesh::err << "ERROR: Bad FEType.radial_family= " << fe_t.radial_family << std::endl; 00361 libmesh_error(); 00362 break; 00363 } 00364 00365 break; 00366 } 00367 00368 00369 00370 00371 // 3D 00372 case 3: 00373 { 00374 switch (fe_t.radial_family) 00375 { 00376 case INFINITE_MAP: 00377 { 00378 libMesh::err << "ERROR: INFINTE_MAP is not a valid shape family for radial approximation." << std::endl; 00379 libmesh_error(); 00380 break; 00381 } 00382 00383 case JACOBI_20_00: 00384 { 00385 switch (fe_t.inf_map) 00386 { 00387 case CARTESIAN: 00388 { 00389 InfFE<3,JACOBI_20_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00390 break; 00391 } 00392 default: 00393 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00394 libmesh_error(); 00395 } 00396 break; 00397 } 00398 00399 case JACOBI_30_00: 00400 { 00401 switch (fe_t.inf_map) 00402 { 00403 case CARTESIAN: 00404 { 00405 InfFE<3,JACOBI_30_00,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00406 break; 00407 } 00408 default: 00409 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00410 libmesh_error(); 00411 } 00412 break; 00413 } 00414 00415 case LEGENDRE: 00416 { 00417 switch (fe_t.inf_map) 00418 { 00419 case CARTESIAN: 00420 { 00421 InfFE<3,LEGENDRE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00422 break; 00423 } 00424 default: 00425 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00426 libmesh_error(); 00427 } 00428 break; 00429 } 00430 00431 case LAGRANGE: 00432 { 00433 switch (fe_t.inf_map) 00434 { 00435 case CARTESIAN: 00436 { 00437 InfFE<3,LAGRANGE,CARTESIAN>::nodal_soln(fe_t, elem, elem_soln, nodal_soln); 00438 break; 00439 } 00440 default: 00441 libMesh::err << "ERROR: Spherical & Ellipsoidal IFEMs not implemented." << std::endl; 00442 libmesh_error(); 00443 } 00444 break; 00445 } 00446 00447 00448 00449 default: 00450 libMesh::err << "ERROR: Bad FEType.radial_family= " << fe_t.radial_family << std::endl; 00451 libmesh_error(); 00452 break; 00453 } 00454 00455 break; 00456 } 00457 00458 default: 00459 libmesh_error(); 00460 } 00461 return; 00462 }
| bool libMesh::FEInterface::ifem_on_reference_element | ( | const Point & | p, | |
| const ElemType | t, | |||
| const Real | eps | |||
| ) | [static, private] |
Definition at line 651 of file fe_interface_inf_fe.C.
References on_reference_element().
00654 { 00655 return FEBase::on_reference_element(p,t,eps); 00656 }
| Real libMesh::FEInterface::ifem_shape | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const Elem * | elem, | |||
| const unsigned int | i, | |||
| const Point & | p | |||
| ) | [static, private] |
Definition at line 766 of file fe_interface_inf_fe.C.
References libMeshEnums::INFINITE_MAP, libMeshEnums::JACOBI_20_00, libMeshEnums::JACOBI_30_00, libMeshEnums::LAGRANGE, libMeshEnums::LEGENDRE, libMesh::FEType::radial_family, and shape().
00771 { 00772 switch (dim) 00773 { 00774 // 1D 00775 case 1: 00776 { 00777 switch (fe_t.radial_family) 00778 { 00779 /* 00780 * For no derivatives (and local coordinates, as 00781 * given in \p p) the infinite element shapes 00782 * are independent of mapping type 00783 */ 00784 case INFINITE_MAP: 00785 return InfFE<1,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p); 00786 00787 case JACOBI_20_00: 00788 return InfFE<1,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p); 00789 00790 case JACOBI_30_00: 00791 return InfFE<1,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p); 00792 00793 case LEGENDRE: 00794 return InfFE<1,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p); 00795 00796 case LAGRANGE: 00797 return InfFE<1,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p); 00798 00799 default: 00800 libmesh_error(); 00801 } 00802 } 00803 00804 00805 // 2D 00806 case 2: 00807 { 00808 switch (fe_t.radial_family) 00809 { 00810 case INFINITE_MAP: 00811 return InfFE<2,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p); 00812 00813 case JACOBI_20_00: 00814 return InfFE<2,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p); 00815 00816 case JACOBI_30_00: 00817 return InfFE<2,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p); 00818 00819 case LEGENDRE: 00820 return InfFE<2,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p); 00821 00822 case LAGRANGE: 00823 return InfFE<2,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p); 00824 00825 default: 00826 libmesh_error(); 00827 } 00828 00829 } 00830 00831 00832 // 3D 00833 case 3: 00834 { 00835 switch (fe_t.radial_family) 00836 { 00837 case INFINITE_MAP: 00838 return InfFE<3,INFINITE_MAP,CARTESIAN>::shape(fe_t, elem, i, p); 00839 00840 case JACOBI_20_00: 00841 return InfFE<3,JACOBI_20_00,CARTESIAN>::shape(fe_t, elem, i, p); 00842 00843 case JACOBI_30_00: 00844 return InfFE<3,JACOBI_30_00,CARTESIAN>::shape(fe_t, elem, i, p); 00845 00846 case LEGENDRE: 00847 return InfFE<3,LEGENDRE,CARTESIAN>::shape(fe_t, elem, i, p); 00848 00849 case LAGRANGE: 00850 return InfFE<3,LAGRANGE,CARTESIAN>::shape(fe_t, elem, i, p); 00851 00852 default: 00853 libmesh_error(); 00854 } 00855 00856 } 00857 00858 00859 default: 00860 libmesh_error(); 00861 } 00862 00863 00864 libmesh_error(); 00865 return 0.; 00866 }
| Real libMesh::FEInterface::ifem_shape | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const ElemType | t, | |||
| const unsigned int | i, | |||
| const Point & | p | |||
| ) | [static, private] |
Definition at line 661 of file fe_interface_inf_fe.C.
References libMeshEnums::INFINITE_MAP, libMeshEnums::JACOBI_20_00, libMeshEnums::JACOBI_30_00, libMeshEnums::LAGRANGE, libMeshEnums::LEGENDRE, libMesh::FEType::radial_family, and shape().
Referenced by shape().
00666 { 00667 switch (dim) 00668 { 00669 // 1D 00670 case 1: 00671 { 00672 switch (fe_t.radial_family) 00673 { 00674 /* 00675 * For no derivatives (and local coordinates, as 00676 * given in \p p) the infinite element shapes 00677 * are independent of mapping type 00678 */ 00679 case INFINITE_MAP: 00680 return InfFE<1,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p); 00681 00682 case JACOBI_20_00: 00683 return InfFE<1,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p); 00684 00685 case JACOBI_30_00: 00686 return InfFE<1,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p); 00687 00688 case LEGENDRE: 00689 return InfFE<1,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p); 00690 00691 case LAGRANGE: 00692 return InfFE<1,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p); 00693 00694 default: 00695 libmesh_error(); 00696 } 00697 } 00698 00699 00700 // 2D 00701 case 2: 00702 { 00703 switch (fe_t.radial_family) 00704 { 00705 case INFINITE_MAP: 00706 return InfFE<2,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p); 00707 00708 case JACOBI_20_00: 00709 return InfFE<2,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p); 00710 00711 case JACOBI_30_00: 00712 return InfFE<2,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p); 00713 00714 case LEGENDRE: 00715 return InfFE<2,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p); 00716 00717 case LAGRANGE: 00718 return InfFE<2,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p); 00719 00720 default: 00721 libmesh_error(); 00722 } 00723 00724 } 00725 00726 00727 // 3D 00728 case 3: 00729 { 00730 switch (fe_t.radial_family) 00731 { 00732 case INFINITE_MAP: 00733 return InfFE<3,INFINITE_MAP,CARTESIAN>::shape(fe_t, t, i, p); 00734 00735 case JACOBI_20_00: 00736 return InfFE<3,JACOBI_20_00,CARTESIAN>::shape(fe_t, t, i, p); 00737 00738 case JACOBI_30_00: 00739 return InfFE<3,JACOBI_30_00,CARTESIAN>::shape(fe_t, t, i, p); 00740 00741 case LEGENDRE: 00742 return InfFE<3,LEGENDRE,CARTESIAN>::shape(fe_t, t, i, p); 00743 00744 case LAGRANGE: 00745 return InfFE<3,LAGRANGE,CARTESIAN>::shape(fe_t, t, i, p); 00746 00747 default: 00748 libmesh_error(); 00749 } 00750 00751 } 00752 00753 00754 default: 00755 libmesh_error(); 00756 } 00757 00758 00759 libmesh_error(); 00760 return 0.; 00761 }
| 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_pointslocated 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 vectorreference_points. The optional parametertolerancedefines how close is "good enough." The map inversion iteration computes the sequence
, and the iteration is terminated when
Definition at line 565 of file fe_interface.C.
References libMesh::err, ifem_inverse_map(), inverse_map(), is_InfFE_elem(), libMesh::TypeVector< T >::size(), and libMesh::Elem::type().
00572 { 00573 const std::size_t n_pts = physical_points.size(); 00574 00575 // Resize the vector 00576 reference_points.resize(n_pts); 00577 00578 if (n_pts == 0) 00579 { 00580 libMesh::err << "WARNING: empty vector physical_points!" 00581 << std::endl; 00582 libmesh_here(); 00583 return; 00584 } 00585 00586 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00587 00588 if ( is_InfFE_elem(elem->type()) ) 00589 { 00590 ifem_inverse_map(dim, fe_t, elem, physical_points, reference_points, tolerance, secure); 00591 return; 00592 00593 // libMesh::err << "ERROR: Not implemented!" 00594 // << std::endl; 00595 // libmesh_error(); 00596 } 00597 00598 #endif 00599 00600 void_fe_with_vec_switch(inverse_map(elem, physical_points, reference_points, tolerance, secure)); 00601 00602 libmesh_error(); 00603 return; 00604 }
| 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
plocated in physical space. This function requires inverting the (probably nonlinear) transformation map, so it is not trivial. The optional parametertolerancedefines how close is "good enough." The map inversion iteration computes the sequence
, and the iteration is terminated when
Definition at line 542 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< OutputType >::coarsened_dof_values(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), libMesh::FEGenericBase< OutputType >::compute_proj_constraints(), libMesh::InfQuad4::contains_point(), libMesh::InfPrism6::contains_point(), libMesh::InfHex8::contains_point(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), ifem_inverse_map(), inverse_map(), 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().
00548 { 00549 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00550 00551 if ( is_InfFE_elem(elem->type()) ) 00552 return ifem_inverse_map(dim, fe_t, elem, p,tolerance, secure); 00553 00554 #endif 00555 00556 fe_with_vec_switch(inverse_map(elem, p, tolerance, secure)); 00557 00558 libmesh_error(); 00559 return Point(); 00560 }
| bool libMesh::FEInterface::is_InfFE_elem | ( | const ElemType | et | ) | [inline, static, private] |
- Returns:
- true if
etis an element to be processed by classInfFE. Otherwise, it returns false. For compatibility with disabled infinite elements it always returns false.
Definition at line 449 of file fe_interface.h.
Referenced by compute_data(), inverse_map(), n_dofs(), n_dofs_at_node(), n_dofs_per_elem(), n_shape_functions(), nodal_soln(), and shape().
| Point libMesh::FEInterface::map | ( | unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const Elem * | elem, | |||
| const Point & | p | |||
| ) | [static] |
Returns the point in physical space of the reference point refpoint which is passed in.
Definition at line 527 of file fe_interface.C.
Referenced by libMesh::Elem::point_test().
00531 { 00532 fe_with_vec_switch(map(elem, p)); 00533 00534 libmesh_error(); 00535 return Point(); 00536 }
| 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 964 of file fe_interface.C.
References libMeshEnums::BERNSTEIN, libMeshEnums::CLOUGH, libMeshEnums::EDGE2, libMeshEnums::EDGE3, libMeshEnums::EDGE4, libMesh::FEType::family, libMeshEnums::HERMITE, libMeshEnums::HEX20, libMeshEnums::HEX27, libMeshEnums::HEX8, libMeshEnums::HIERARCHIC, libMeshEnums::L2_HIERARCHIC, libMeshEnums::L2_LAGRANGE, libMeshEnums::LAGRANGE, libMeshEnums::LAGRANGE_VEC, libMeshEnums::MONOMIAL, libMeshEnums::NEDELEC_ONE, libMeshEnums::PRISM15, libMeshEnums::PRISM18, libMeshEnums::PRISM6, libMeshEnums::PYRAMID5, libMeshEnums::QUAD4, libMeshEnums::QUAD8, libMeshEnums::QUAD9, libMeshEnums::SZABAB, libMeshEnums::TET10, libMeshEnums::TET4, libMeshEnums::TRI3, libMeshEnums::TRI6, and libMeshEnums::XYZ.
Referenced by libMesh::DofMap::reinit().
00966 { 00967 // Yeah, I know, infinity is much larger than 11, but our 00968 // solvers don't seem to like high degree polynomials, and our 00969 // quadrature rules and number_lookups tables 00970 // need to go up higher. 00971 const unsigned int unlimited = 11; 00972 00973 // If we used 0 as a default, then elements missing from this 00974 // table (e.g. infinite elements) would be considered broken. 00975 const unsigned int unknown = unlimited; 00976 00977 switch (fe_t.family) 00978 { 00979 case LAGRANGE: 00980 case L2_LAGRANGE: // TODO: L2_LAGRANGE can have higher "max_order" than LAGRANGE 00981 case LAGRANGE_VEC: 00982 switch (el_t) 00983 { 00984 case EDGE2: 00985 case EDGE3: 00986 case EDGE4: 00987 return 3; 00988 case TRI3: 00989 return 1; 00990 case TRI6: 00991 return 2; 00992 case QUAD4: 00993 return 1; 00994 case QUAD8: 00995 case QUAD9: 00996 return 2; 00997 case TET4: 00998 return 1; 00999 case TET10: 01000 return 2; 01001 case HEX8: 01002 return 1; 01003 case HEX20: 01004 case HEX27: 01005 return 2; 01006 case PRISM6: 01007 case PRISM15: 01008 return 1; 01009 case PRISM18: 01010 return 2; 01011 case PYRAMID5: 01012 return 1; 01013 default: 01014 return unknown; 01015 } 01016 break; 01017 case MONOMIAL: 01018 switch (el_t) 01019 { 01020 case EDGE2: 01021 case EDGE3: 01022 case EDGE4: 01023 case TRI3: 01024 case TRI6: 01025 case QUAD4: 01026 case QUAD8: 01027 case QUAD9: 01028 case TET4: 01029 case TET10: 01030 case HEX8: 01031 case HEX20: 01032 case HEX27: 01033 case PRISM6: 01034 case PRISM15: 01035 case PRISM18: 01036 case PYRAMID5: 01037 return unlimited; 01038 default: 01039 return unknown; 01040 } 01041 break; 01042 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 01043 case BERNSTEIN: 01044 switch (el_t) 01045 { 01046 case EDGE2: 01047 case EDGE3: 01048 case EDGE4: 01049 return unlimited; 01050 case TRI3: 01051 return 0; 01052 case TRI6: 01053 return 6; 01054 case QUAD4: 01055 return 0; 01056 case QUAD8: 01057 case QUAD9: 01058 return unlimited; 01059 case TET4: 01060 return 1; 01061 case TET10: 01062 return 2; 01063 case HEX8: 01064 return 0; 01065 case HEX20: 01066 return 2; 01067 case HEX27: 01068 return 4; 01069 case PRISM6: 01070 case PRISM15: 01071 case PRISM18: 01072 case PYRAMID5: 01073 return 0; 01074 default: 01075 return unknown; 01076 } 01077 break; 01078 case SZABAB: 01079 switch (el_t) 01080 { 01081 case EDGE2: 01082 case EDGE3: 01083 case EDGE4: 01084 return 7; 01085 case TRI3: 01086 return 0; 01087 case TRI6: 01088 return 7; 01089 case QUAD4: 01090 return 0; 01091 case QUAD8: 01092 case QUAD9: 01093 return 7; 01094 case TET4: 01095 case TET10: 01096 case HEX8: 01097 case HEX20: 01098 case HEX27: 01099 case PRISM6: 01100 case PRISM15: 01101 case PRISM18: 01102 case PYRAMID5: 01103 return 0; 01104 default: 01105 return unknown; 01106 } 01107 break; 01108 #endif 01109 case XYZ: 01110 switch (el_t) 01111 { 01112 case EDGE2: 01113 case EDGE3: 01114 case EDGE4: 01115 case TRI3: 01116 case TRI6: 01117 case QUAD4: 01118 case QUAD8: 01119 case QUAD9: 01120 case TET4: 01121 case TET10: 01122 case HEX8: 01123 case HEX20: 01124 case HEX27: 01125 case PRISM6: 01126 case PRISM15: 01127 case PRISM18: 01128 case PYRAMID5: 01129 return unlimited; 01130 default: 01131 return unknown; 01132 } 01133 break; 01134 case CLOUGH: 01135 switch (el_t) 01136 { 01137 case EDGE2: 01138 case EDGE3: 01139 return 3; 01140 case EDGE4: 01141 case TRI3: 01142 return 0; 01143 case TRI6: 01144 return 3; 01145 case QUAD4: 01146 case QUAD8: 01147 case QUAD9: 01148 case TET4: 01149 case TET10: 01150 case HEX8: 01151 case HEX20: 01152 case HEX27: 01153 case PRISM6: 01154 case PRISM15: 01155 case PRISM18: 01156 case PYRAMID5: 01157 return 0; 01158 default: 01159 return unknown; 01160 } 01161 break; 01162 case HERMITE: 01163 switch (el_t) 01164 { 01165 case EDGE2: 01166 case EDGE3: 01167 return unlimited; 01168 case EDGE4: 01169 case TRI3: 01170 case TRI6: 01171 return 0; 01172 case QUAD4: 01173 return 3; 01174 case QUAD8: 01175 case QUAD9: 01176 return unlimited; 01177 case TET4: 01178 case TET10: 01179 return 0; 01180 case HEX8: 01181 return 3; 01182 case HEX20: 01183 case HEX27: 01184 return unlimited; 01185 case PRISM6: 01186 case PRISM15: 01187 case PRISM18: 01188 case PYRAMID5: 01189 return 0; 01190 default: 01191 return unknown; 01192 } 01193 break; 01194 case HIERARCHIC: 01195 switch (el_t) 01196 { 01197 case EDGE2: 01198 case EDGE3: 01199 case EDGE4: 01200 return unlimited; 01201 case TRI3: 01202 return 1; 01203 case TRI6: 01204 return unlimited; 01205 case QUAD4: 01206 return 1; 01207 case QUAD8: 01208 case QUAD9: 01209 return unlimited; 01210 case TET4: 01211 case TET10: 01212 return 0; 01213 case HEX8: 01214 case HEX20: 01215 return 1; 01216 case HEX27: 01217 return unlimited; 01218 case PRISM6: 01219 case PRISM15: 01220 case PRISM18: 01221 case PYRAMID5: 01222 return 0; 01223 default: 01224 return unknown; 01225 } 01226 break; 01227 case L2_HIERARCHIC: 01228 switch (el_t) 01229 { 01230 case EDGE2: 01231 case EDGE3: 01232 case EDGE4: 01233 return unlimited; 01234 case TRI3: 01235 return 1; 01236 case TRI6: 01237 return unlimited; 01238 case QUAD4: 01239 return 1; 01240 case QUAD8: 01241 case QUAD9: 01242 return unlimited; 01243 case TET4: 01244 case TET10: 01245 return 0; 01246 case HEX8: 01247 case HEX20: 01248 return 1; 01249 case HEX27: 01250 return unlimited; 01251 case PRISM6: 01252 case PRISM15: 01253 case PRISM18: 01254 case PYRAMID5: 01255 return 0; 01256 default: 01257 return unknown; 01258 } 01259 break; 01260 case NEDELEC_ONE: 01261 switch (el_t) 01262 { 01263 case TRI6: 01264 case QUAD8: 01265 case QUAD9: 01266 return 1; 01267 default: 01268 return 0; 01269 } 01270 break; 01271 default: 01272 return 0; 01273 break; 01274 } 01275 }
| 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 404 of file fe_interface.C.
References ifem_n_dofs(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), compute_data(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), ifem_n_dofs(), and libMesh::HPCoarsenTest::select_refinement().
00407 { 00408 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00409 00410 if ( is_InfFE_elem(t) ) 00411 return ifem_n_dofs(dim, fe_t, t); 00412 00413 #endif 00414 00415 const Order o = fe_t.order; 00416 00417 fe_with_vec_switch(n_dofs(t, o)); 00418 00419 libmesh_error(); 00420 return 0; 00421 }
| 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 426 of file fe_interface.C.
References ifem_n_dofs_at_node(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::DofMap::constrain_p_dofs(), libMesh::DofMap::dof_indices(), ifem_n_dofs_at_node(), libMesh::DofMap::old_dof_indices(), libMesh::BoundaryProjectSolution::operator()(), libMesh::ProjectFEMSolution::operator()(), libMesh::ProjectSolution::operator()(), and libMesh::DofMap::reinit().
00430 { 00431 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00432 00433 if ( is_InfFE_elem(t) ) 00434 return ifem_n_dofs_at_node(dim, fe_t, t, n); 00435 00436 #endif 00437 00438 const Order o = fe_t.order; 00439 00440 fe_with_vec_switch(n_dofs_at_node(t, o, n)); 00441 00442 libmesh_error(); 00443 return 0; 00444 }
| 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 450 of file fe_interface.C.
References ifem_n_dofs_per_elem(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by libMesh::DofMap::dof_indices(), ifem_n_dofs_per_elem(), libMesh::DofMap::old_dof_indices(), and libMesh::DofMap::reinit().
00453 { 00454 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00455 00456 if ( is_InfFE_elem(t) ) 00457 return ifem_n_dofs_per_elem(dim, fe_t, t); 00458 00459 #endif 00460 00461 const Order o = fe_t.order; 00462 00463 fe_with_vec_switch(n_dofs_per_elem(t, o)); 00464 00465 libmesh_error(); 00466 return 0; 00467 }
| 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 375 of file fe_interface.C.
References ifem_n_shape_functions(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by ifem_n_shape_functions().
00378 { 00379 00380 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00381 /* 00382 * Since the FEType, stored in DofMap/(some System child), has to 00383 * be the _same_ for InfFE and FE, we have to catch calls 00384 * to infinite elements through the element type. 00385 */ 00386 00387 if ( is_InfFE_elem(t) ) 00388 return ifem_n_shape_functions(dim, fe_t, t); 00389 00390 #endif 00391 00392 const Order o = fe_t.order; 00393 00394 fe_with_vec_switch(n_shape_functions(t, o)); 00395 00396 libmesh_error(); 00397 return 0; 00398 }
| 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 1320 of file fe_interface.C.
References libMesh::FEType::family, libMeshEnums::LAGRANGE_VEC, libMesh::MeshBase::mesh_dimension(), and libMeshEnums::NEDELEC_ONE.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::EquationSystems::build_solution_vector(), and libMesh::EquationSystems::build_variable_names().
01321 { 01322 switch (fe_type.family) 01323 { 01324 //FIXME: We currently assume that the number of vector components is tied 01325 // to the mesh dimension. This will break for mixed-dimension meshes. 01326 case LAGRANGE_VEC: 01327 case NEDELEC_ONE: 01328 return mesh.mesh_dimension(); 01329 default: 01330 return 1; 01331 } 01332 }
| 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 503 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_solution_vector(), ifem_nodal_soln(), libMesh::EnsightIO::write_scalar_ascii(), and libMesh::EnsightIO::write_vector_ascii().
00508 { 00509 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00510 00511 if ( is_InfFE_elem(elem->type()) ) 00512 { 00513 ifem_nodal_soln(dim, fe_t, elem, elem_soln, nodal_soln); 00514 return; 00515 } 00516 00517 #endif 00518 00519 const Order order = fe_t.order; 00520 00521 void_fe_with_vec_switch(nodal_soln(elem, order, elem_soln, nodal_soln)); 00522 }
| 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 here the parameter eps can be specified to indicate a tolerance. For example,
becomes
.
Definition at line 608 of file fe_interface.C.
Referenced by libMesh::InfQuad4::contains_point(), libMesh::InfPrism6::contains_point(), libMesh::InfHex8::contains_point(), ifem_on_reference_element(), and libMesh::Elem::point_test().
00611 { 00612 return FEBase::on_reference_element(p,t,eps); 00613 }
| void libMesh::FEInterface::shape | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const Elem * | elem, | |||
| const unsigned int | i, | |||
| const Point & | p, | |||
| RealGradient & | phi | |||
| ) | [inline] |
| void libMesh::FEInterface::shape | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const ElemType | t, | |||
| const unsigned int | i, | |||
| const Point & | p, | |||
| RealGradient & | phi | |||
| ) | [inline] |
| void libMesh::FEInterface::shape | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const Elem * | elem, | |||
| const unsigned int | i, | |||
| const Point & | p, | |||
| Real & | phi | |||
| ) | [inline] |
| void libMesh::FEInterface::shape | ( | const unsigned int | dim, | |
| const FEType & | fe_t, | |||
| const ElemType | t, | |||
| const unsigned int | i, | |||
| const Point & | p, | |||
| Real & | phi | |||
| ) | [inline] |
| 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 | |||
| ) | [inline, static] |
- Returns:
- the value of the
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.
On a p-refined element, fe_t.order should be the total order of the element.
| 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 | |||
| ) | [inline, static] |
- Returns:
- the value of the
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.
On a p-refined element, fe_t.order should be the total order of the element.
| 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
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.
On a p-refined element, fe_t.order should be the base order of the element.
Definition at line 639 of file fe_interface.C.
References ifem_shape(), is_InfFE_elem(), libMesh::FEType::order, shape(), and libMesh::Elem::type().
00644 { 00645 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00646 00647 if ( is_InfFE_elem(elem->type()) ) 00648 return ifem_shape(dim, fe_t, elem, i, p); 00649 00650 #endif 00651 00652 const Order o = fe_t.order; 00653 00654 fe_switch(shape(elem,o,i,p)); 00655 00656 libmesh_error(); 00657 return 0.; 00658 }
| 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
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.
On a p-refined element, fe_t.order should be the total order of the element.
Definition at line 618 of file fe_interface.C.
References ifem_shape(), is_InfFE_elem(), and libMesh::FEType::order.
Referenced by compute_data(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), ifem_shape(), and shape().
00623 { 00624 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00625 00626 if ( is_InfFE_elem(t) ) 00627 return ifem_shape(dim, fe_t, t, i, p); 00628 00629 #endif 00630 00631 const Order o = fe_t.order; 00632 00633 fe_switch(shape(t,o,i,p)); 00634 00635 libmesh_error(); 00636 return 0.; 00637 }
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:24 UTC
Hosted By: