libMesh
|
Class contained in FE that encapsulates mapping (i.e. More...
#include <fe_map.h>
Public Member Functions | |
FEMap (Real jtol=0) | |
virtual | ~FEMap ()=default |
template<unsigned int Dim> | |
void | init_reference_to_physical_map (const std::vector< Point > &qp, const Elem *elem) |
void | init_reference_to_physical_map (unsigned int dim, const std::vector< Point > &qp, const Elem *elem) |
void | compute_single_point_map (const unsigned int dim, const std::vector< Real > &qw, const Elem *elem, unsigned int p, const std::vector< const Node *> &elem_nodes, bool compute_second_derivatives) |
Compute the jacobian and some other additional data fields at the single point with index p. More... | |
virtual void | compute_affine_map (const unsigned int dim, const std::vector< Real > &qw, const Elem *elem) |
Compute the jacobian and some other additional data fields. More... | |
virtual void | compute_null_map (const unsigned int dim, const std::vector< Real > &qw) |
Assign a fake jacobian and some other additional data fields. More... | |
virtual void | compute_map (const unsigned int dim, const std::vector< Real > &qw, const Elem *elem, bool calculate_d2phi) |
Compute the jacobian and some other additional data fields. More... | |
virtual void | compute_face_map (int dim, const std::vector< Real > &qw, const Elem *side) |
Same as compute_map, but for a side. More... | |
void | compute_edge_map (int dim, const std::vector< Real > &qw, const Elem *side) |
Same as before, but for an edge. More... | |
template<unsigned int Dim> | |
void | init_face_shape_functions (const std::vector< Point > &qp, const Elem *side) |
Initializes the reference to physical element map for a side. More... | |
template<unsigned int Dim> | |
void | init_edge_shape_functions (const std::vector< Point > &qp, const Elem *edge) |
Same as before, but for an edge. More... | |
const std::vector< Point > & | get_xyz () const |
const std::vector< Real > & | get_jacobian () const |
const std::vector< Real > & | get_JxW () const |
const std::vector< RealGradient > & | get_dxyzdxi () const |
const std::vector< RealGradient > & | get_dxyzdeta () const |
const std::vector< RealGradient > & | get_dxyzdzeta () const |
const std::vector< RealGradient > & | get_d2xyzdxi2 () const |
const std::vector< RealGradient > & | get_d2xyzdeta2 () const |
const std::vector< RealGradient > & | get_d2xyzdzeta2 () const |
const std::vector< RealGradient > & | get_d2xyzdxideta () const |
const std::vector< RealGradient > & | get_d2xyzdxidzeta () const |
const std::vector< RealGradient > & | get_d2xyzdetadzeta () const |
const std::vector< Real > & | get_dxidx () const |
const std::vector< Real > & | get_dxidy () const |
const std::vector< Real > & | get_dxidz () const |
const std::vector< Real > & | get_detadx () const |
const std::vector< Real > & | get_detady () const |
const std::vector< Real > & | get_detadz () const |
const std::vector< Real > & | get_dzetadx () const |
const std::vector< Real > & | get_dzetady () const |
const std::vector< Real > & | get_dzetadz () const |
const std::vector< std::vector< Real > > & | get_d2xidxyz2 () const |
Second derivatives of "xi" reference coordinate wrt physical coordinates. More... | |
const std::vector< std::vector< Real > > & | get_d2etadxyz2 () const |
Second derivatives of "eta" reference coordinate wrt physical coordinates. More... | |
const std::vector< std::vector< Real > > & | get_d2zetadxyz2 () const |
Second derivatives of "zeta" reference coordinate wrt physical coordinates. More... | |
const std::vector< std::vector< Real > > & | get_psi () const |
const std::vector< std::vector< Real > > & | get_phi_map () const |
const std::vector< std::vector< Real > > & | get_dphidxi_map () const |
const std::vector< std::vector< Real > > & | get_dphideta_map () const |
const std::vector< std::vector< Real > > & | get_dphidzeta_map () const |
const std::vector< std::vector< Point > > & | get_tangents () const |
const std::vector< Point > & | get_normals () const |
const std::vector< Real > & | get_curvatures () const |
void | print_JxW (std::ostream &os) const |
Prints the Jacobian times the weight for each quadrature point. More... | |
void | print_xyz (std::ostream &os) const |
Prints the spatial location of each quadrature point (on the physical element). More... | |
std::vector< std::vector< Real > > & | get_psi () |
std::vector< std::vector< Real > > & | get_dpsidxi () |
const std::vector< std::vector< Real > > & | get_dpsidxi () const |
std::vector< std::vector< Real > > & | get_dpsideta () |
const std::vector< std::vector< Real > > & | get_dpsideta () const |
std::vector< std::vector< Real > > & | get_d2psidxi2 () |
const std::vector< std::vector< Real > > & | get_d2psidxi2 () const |
std::vector< std::vector< Real > > & | get_d2psidxideta () |
const std::vector< std::vector< Real > > & | get_d2psidxideta () const |
std::vector< std::vector< Real > > & | get_d2psideta2 () |
const std::vector< std::vector< Real > > & | get_d2psideta2 () const |
std::vector< std::vector< Real > > & | get_phi_map () |
std::vector< std::vector< Real > > & | get_dphidxi_map () |
std::vector< std::vector< Real > > & | get_dphideta_map () |
std::vector< std::vector< Real > > & | get_dphidzeta_map () |
std::vector< std::vector< Real > > & | get_d2phidxi2_map () |
std::vector< std::vector< Real > > & | get_d2phidxideta_map () |
std::vector< std::vector< Real > > & | get_d2phidxidzeta_map () |
std::vector< std::vector< Real > > & | get_d2phideta2_map () |
std::vector< std::vector< Real > > & | get_d2phidetadzeta_map () |
std::vector< std::vector< Real > > & | get_d2phidzeta2_map () |
std::vector< Real > & | get_JxW () |
void | add_calculations () |
Allows the user to prerequest additional calculations in between two calls to reinit();. More... | |
void | set_jacobian_tolerance (Real tol) |
Set the Jacobian tolerance used for determining when the mapping fails. More... | |
Static Public Member Functions | |
static std::unique_ptr< FEMap > | build (FEType fe_type) |
static FEFamily | map_fe_type (const Elem &elem) |
static Point | map (const unsigned int dim, const Elem *elem, const Point &reference_point) |
static Point | map_deriv (const unsigned int dim, const Elem *elem, const unsigned int j, const Point &reference_point) |
static Point | inverse_map (const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true) |
static void | inverse_map (unsigned int dim, const Elem *elem, const std::vector< Point > &physical_points, std::vector< Point > &reference_points, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true) |
Takes a number points in physical space (in the physical_points vector) and finds their location on the reference element for the input element elem . More... | |
Protected Member Functions | |
void | determine_calculations () |
Determine which values are to be calculated. More... | |
void | resize_quadrature_map_vectors (const unsigned int dim, unsigned int n_qp) |
A utility function for use by compute_*_map. More... | |
Real | dxdxi_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Real | dydxi_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Real | dzdxi_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Real | dxdeta_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Real | dydeta_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Real | dzdeta_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Real | dxdzeta_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Real | dydzeta_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Real | dzdzeta_map (const unsigned int p) const |
Used in FEMap::compute_map() , which should be be usable in derived classes, and therefore protected. More... | |
Protected Attributes | |
std::vector< Point > | xyz |
The spatial locations of the quadrature points. More... | |
std::vector< RealGradient > | dxyzdxi_map |
Vector of partial derivatives: d(x)/d(xi), d(y)/d(xi), d(z)/d(xi) More... | |
std::vector< RealGradient > | dxyzdeta_map |
Vector of partial derivatives: d(x)/d(eta), d(y)/d(eta), d(z)/d(eta) More... | |
std::vector< RealGradient > | dxyzdzeta_map |
Vector of partial derivatives: d(x)/d(zeta), d(y)/d(zeta), d(z)/d(zeta) More... | |
std::vector< RealGradient > | d2xyzdxi2_map |
Vector of second partial derivatives in xi: d^2(x)/d(xi)^2, d^2(y)/d(xi)^2, d^2(z)/d(xi)^2. More... | |
std::vector< RealGradient > | d2xyzdxideta_map |
Vector of mixed second partial derivatives in xi-eta: d^2(x)/d(xi)d(eta) d^2(y)/d(xi)d(eta) d^2(z)/d(xi)d(eta) More... | |
std::vector< RealGradient > | d2xyzdeta2_map |
Vector of second partial derivatives in eta: d^2(x)/d(eta)^2. More... | |
std::vector< RealGradient > | d2xyzdxidzeta_map |
Vector of second partial derivatives in xi-zeta: d^2(x)/d(xi)d(zeta), d^2(y)/d(xi)d(zeta), d^2(z)/d(xi)d(zeta) More... | |
std::vector< RealGradient > | d2xyzdetadzeta_map |
Vector of mixed second partial derivatives in eta-zeta: d^2(x)/d(eta)d(zeta) d^2(y)/d(eta)d(zeta) d^2(z)/d(eta)d(zeta) More... | |
std::vector< RealGradient > | d2xyzdzeta2_map |
Vector of second partial derivatives in zeta: d^2(x)/d(zeta)^2. More... | |
std::vector< Real > | dxidx_map |
Map for partial derivatives: d(xi)/d(x). More... | |
std::vector< Real > | dxidy_map |
Map for partial derivatives: d(xi)/d(y). More... | |
std::vector< Real > | dxidz_map |
Map for partial derivatives: d(xi)/d(z). More... | |
std::vector< Real > | detadx_map |
Map for partial derivatives: d(eta)/d(x). More... | |
std::vector< Real > | detady_map |
Map for partial derivatives: d(eta)/d(y). More... | |
std::vector< Real > | detadz_map |
Map for partial derivatives: d(eta)/d(z). More... | |
std::vector< Real > | dzetadx_map |
Map for partial derivatives: d(zeta)/d(x). More... | |
std::vector< Real > | dzetady_map |
Map for partial derivatives: d(zeta)/d(y). More... | |
std::vector< Real > | dzetadz_map |
Map for partial derivatives: d(zeta)/d(z). More... | |
std::vector< std::vector< Real > > | d2xidxyz2_map |
Second derivatives of "xi" reference coordinate wrt physical coordinates. More... | |
std::vector< std::vector< Real > > | d2etadxyz2_map |
Second derivatives of "eta" reference coordinate wrt physical coordinates. More... | |
std::vector< std::vector< Real > > | d2zetadxyz2_map |
Second derivatives of "zeta" reference coordinate wrt physical coordinates. More... | |
std::vector< std::vector< Real > > | phi_map |
Map for the shape function phi. More... | |
std::vector< std::vector< Real > > | dphidxi_map |
Map for the derivative, d(phi)/d(xi). More... | |
std::vector< std::vector< Real > > | dphideta_map |
Map for the derivative, d(phi)/d(eta). More... | |
std::vector< std::vector< Real > > | dphidzeta_map |
Map for the derivative, d(phi)/d(zeta). More... | |
std::vector< std::vector< Real > > | d2phidxi2_map |
Map for the second derivative, d^2(phi)/d(xi)^2. More... | |
std::vector< std::vector< Real > > | d2phidxideta_map |
Map for the second derivative, d^2(phi)/d(xi)d(eta). More... | |
std::vector< std::vector< Real > > | d2phidxidzeta_map |
Map for the second derivative, d^2(phi)/d(xi)d(zeta). More... | |
std::vector< std::vector< Real > > | d2phideta2_map |
Map for the second derivative, d^2(phi)/d(eta)^2. More... | |
std::vector< std::vector< Real > > | d2phidetadzeta_map |
Map for the second derivative, d^2(phi)/d(eta)d(zeta). More... | |
std::vector< std::vector< Real > > | d2phidzeta2_map |
Map for the second derivative, d^2(phi)/d(zeta)^2. More... | |
std::vector< std::vector< Real > > | psi_map |
Map for the side shape functions, psi. More... | |
std::vector< std::vector< Real > > | dpsidxi_map |
Map for the derivative of the side functions, d(psi)/d(xi). More... | |
std::vector< std::vector< Real > > | dpsideta_map |
Map for the derivative of the side function, d(psi)/d(eta). More... | |
std::vector< std::vector< Real > > | d2psidxi2_map |
Map for the second derivatives (in xi) of the side shape functions. More... | |
std::vector< std::vector< Real > > | d2psidxideta_map |
Map for the second (cross) derivatives in xi, eta of the side shape functions. More... | |
std::vector< std::vector< Real > > | d2psideta2_map |
Map for the second derivatives (in eta) of the side shape functions. More... | |
std::vector< std::vector< Point > > | tangents |
Tangent vectors on boundary at quadrature points. More... | |
std::vector< Point > | normals |
Normal vectors on boundary at quadrature points. More... | |
std::vector< Real > | curvatures |
The mean curvature (= one half the sum of the principal curvatures) on the boundary at the quadrature points. More... | |
std::vector< Real > | jac |
Jacobian values at quadrature points. More... | |
std::vector< Real > | JxW |
Jacobian*Weight values at quadrature points. More... | |
bool | calculations_started |
Have calculations with this object already been started? Then all get_* functions should already have been called. More... | |
bool | calculate_xyz |
Should we calculate physical point locations? More... | |
bool | calculate_dxyz |
Should we calculate mapping gradients? More... | |
bool | calculate_d2xyz |
Should we calculate mapping hessians? More... | |
Real | jacobian_tolerance |
The Jacobian tolerance used for determining when the mapping fails. More... | |
Private Member Functions | |
void | compute_inverse_map_second_derivs (unsigned p) |
A helper function used by FEMap::compute_single_point_map() to compute second derivatives of the inverse map. More... | |
Private Attributes | |
std::vector< const Node * > | _elem_nodes |
Work vector for compute_affine_map() More... | |
Class contained in FE that encapsulates mapping (i.e.
from physical space to reference space and vice-versa) quantities and operations.
libMesh::FEMap::FEMap | ( | Real | jtol = 0 | ) |
Definition at line 62 of file fe_map.C.
|
virtualdefault |
void libMesh::FEMap::add_calculations | ( | ) |
Allows the user to prerequest additional calculations in between two calls to reinit();.
Definition at line 88 of file fe_map.C.
References calculations_started, d2phideta2_map, d2phidetadzeta_map, d2phidxi2_map, d2phidxideta_map, d2phidxidzeta_map, d2phidzeta2_map, dphideta_map, dphidxi_map, dphidzeta_map, and phi_map.
Definition at line 74 of file fe_map.C.
References libMesh::FEType::family, and libMesh::XYZ.
|
virtual |
Compute the jacobian and some other additional data fields.
Takes the integration weights as input, along with a pointer to the element. The element is assumed to have a constant Jacobian
Definition at line 1281 of file fe_map.C.
References _elem_nodes, calculate_d2xyz, calculate_dxyz, calculate_xyz, compute_single_point_map(), d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, detadx_map, detady_map, detadz_map, dim, dxidx_map, dxidy_map, dxidz_map, dxyzdeta_map, dxyzdxi_map, dxyzdzeta_map, dzetadx_map, dzetady_map, dzetadz_map, libMesh::index_range(), jac, JxW, libMesh::libmesh_assert(), n_nodes, libMesh::Elem::n_nodes(), libMesh::Elem::node_ptr(), phi_map, resize_quadrature_map_vectors(), and xyz.
Referenced by compute_map().
void libMesh::FEMap::compute_edge_map | ( | int | dim, |
const std::vector< Real > & | qw, | ||
const Elem * | side | ||
) |
Same as before, but for an edge.
Useful for some projections.
Definition at line 1005 of file fe_boundary.C.
References calculate_d2xyz, calculate_dxyz, calculate_xyz, compute_face_map(), curvatures, d2psidxi2_map, d2xyzdeta2_map, d2xyzdxi2_map, d2xyzdxideta_map, determine_calculations(), dim, dpsidxi_map, dxdxi_map(), dxyzdeta_map, dxyzdxi_map, dydxi_map(), dzdxi_map(), JxW, libMesh::libmesh_assert(), normals, libMesh::Elem::point(), psi_map, libMesh::Real, std::sqrt(), tangents, and xyz.
|
virtual |
Same as compute_map, but for a side.
Useful for boundary integration.
Reimplemented in libMesh::FEXYZMap.
Definition at line 641 of file fe_boundary.C.
References calculate_d2xyz, calculate_dxyz, calculate_xyz, libMesh::TypeVector< T >::cross(), curvatures, d2psideta2_map, d2psidxi2_map, d2psidxideta_map, d2xyzdeta2_map, d2xyzdxi2_map, d2xyzdxideta_map, libMesh::Elem::default_order(), determine_calculations(), dim, dpsideta_map, dpsidxi_map, dxdeta_map(), dxdxi_map(), dxyzdeta_map, dxyzdxi_map, dydeta_map(), dydxi_map(), dzdeta_map(), dzdxi_map(), libMesh::Elem::interior_parent(), inverse_map(), JxW, libMesh::libmesh_assert(), map_deriv(), map_fe_type(), libMesh::FEInterface::n_shape_functions(), libMesh::Elem::node_id(), normals, libMesh::Elem::point(), psi_map, libMesh::Real, std::sqrt(), tangents, libMesh::TypeVector< T >::unit(), and xyz.
Referenced by compute_edge_map(), and libMesh::FEXYZMap::compute_face_map().
|
private |
A helper function used by FEMap::compute_single_point_map() to compute second derivatives of the inverse map.
Definition at line 1506 of file fe_map.C.
References d2etadxyz2_map, d2xidxyz2_map, d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, d2zetadxyz2_map, detadx_map, detady_map, detadz_map, dxidx_map, dxidy_map, dxidz_map, dzetadx_map, dzetady_map, dzetadz_map, and libMesh::libmesh_ignore().
Referenced by compute_single_point_map().
|
virtual |
Compute the jacobian and some other additional data fields.
Takes the integration weights as input, along with a pointer to the element. Also takes a boolean parameter indicating whether second derivatives need to be calculated, allowing us to potentially skip unnecessary, expensive computations.
Definition at line 1437 of file fe_map.C.
References _elem_nodes, compute_affine_map(), compute_null_map(), compute_single_point_map(), dim, libMesh::MeshTools::Subdivision::find_one_ring(), libMesh::Elem::has_affine_map(), libMesh::libmesh_assert(), libMesh::Elem::n_nodes(), libMesh::Elem::node_index_range(), libMesh::Elem::node_ptr(), resize_quadrature_map_vectors(), libMesh::TRI3SUBDIVISION, and libMesh::Elem::type().
Referenced by libMesh::Elem::has_invertible_map().
|
virtual |
Assign a fake jacobian and some other additional data fields.
Takes the integration weights as input. For use on non-element evaluations.
Definition at line 1362 of file fe_map.C.
References calculate_d2xyz, calculate_dxyz, calculate_xyz, d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, detadx_map, detady_map, detadz_map, dim, dxidx_map, dxidy_map, dxidz_map, dxyzdeta_map, dxyzdxi_map, dxyzdzeta_map, dzetadx_map, dzetady_map, dzetadz_map, jac, JxW, resize_quadrature_map_vectors(), and xyz.
Referenced by compute_map().
void libMesh::FEMap::compute_single_point_map | ( | const unsigned int | dim, |
const std::vector< Real > & | qw, | ||
const Elem * | elem, | ||
unsigned int | p, | ||
const std::vector< const Node *> & | elem_nodes, | ||
bool | compute_second_derivatives | ||
) |
Compute the jacobian and some other additional data fields at the single point with index p.
Takes the integration weights as input, along with a pointer to the element and a list of points that contribute to the element. Also takes a boolean flag telling whether second derivatives should actually be computed.
Definition at line 447 of file fe_map.C.
References calculate_d2xyz, calculate_dxyz, calculate_xyz, calculations_started, compute_inverse_map_second_derivs(), d2etadxyz2_map, d2phideta2_map, d2phidetadzeta_map, d2phidxi2_map, d2phidxideta_map, d2phidxidzeta_map, d2phidzeta2_map, d2xidxyz2_map, d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, d2zetadxyz2_map, detadx_map, detady_map, detadz_map, dim, dphideta_map, dphidxi_map, dphidzeta_map, dxdeta_map(), dxdxi_map(), dxdzeta_map(), dxidx_map, dxidy_map, dxidz_map, dxyzdeta_map, dxyzdxi_map, dxyzdzeta_map, dydeta_map(), dydxi_map(), dydzeta_map(), dzdeta_map(), dzdxi_map(), dzdzeta_map(), dzetadx_map, dzetady_map, dzetadz_map, libMesh::err, libMesh::DofObject::id(), libMesh::index_range(), jac, jacobian_tolerance, JxW, libMesh::libmesh_assert(), libMesh::libmesh_ignore(), phi_map, libMesh::Elem::print_info(), libMesh::Real, std::sqrt(), libMesh::DenseMatrix< T >::vector_mult(), and xyz.
Referenced by compute_affine_map(), and compute_map().
|
inlineprotected |
Determine which values are to be calculated.
Definition at line 648 of file fe_map.h.
References calculate_d2xyz, calculate_dxyz, and calculations_started.
Referenced by compute_edge_map(), compute_face_map(), init_edge_shape_functions(), init_face_shape_functions(), init_reference_to_physical_map(), and resize_quadrature_map_vectors().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 695 of file fe_map.h.
References dxyzdeta_map.
Referenced by compute_face_map(), and compute_single_point_map().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 671 of file fe_map.h.
References dxyzdxi_map.
Referenced by compute_edge_map(), compute_face_map(), and compute_single_point_map().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 719 of file fe_map.h.
References dxyzdzeta_map.
Referenced by compute_single_point_map().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 703 of file fe_map.h.
References dxyzdeta_map.
Referenced by compute_face_map(), and compute_single_point_map().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 679 of file fe_map.h.
References dxyzdxi_map.
Referenced by compute_edge_map(), compute_face_map(), and compute_single_point_map().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 727 of file fe_map.h.
References dxyzdzeta_map.
Referenced by compute_single_point_map().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 711 of file fe_map.h.
References dxyzdeta_map.
Referenced by compute_face_map(), and compute_single_point_map().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 687 of file fe_map.h.
References dxyzdxi_map.
Referenced by compute_edge_map(), compute_face_map(), and compute_single_point_map().
Used in FEMap::compute_map()
, which should be be usable in derived classes, and therefore protected.
Definition at line 735 of file fe_map.h.
References dxyzdzeta_map.
Referenced by compute_single_point_map().
|
inline |
Definition at line 453 of file fe_map.h.
References calculate_d2xyz, calculations_started, curvatures, and libMesh::libmesh_assert().
|
inline |
Second derivatives of "eta" reference coordinate wrt physical coordinates.
Definition at line 388 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2etadxyz2_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_d2phi().
|
inline |
Definition at line 602 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2phideta2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 609 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2phidetadzeta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 581 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2phidxi2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 588 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2phidxideta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 595 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2phidxidzeta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 616 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2phidzeta2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 535 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2psideta2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 543 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2psideta2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 507 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2psidxi2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 514 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2psidxi2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 521 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2psidxideta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 528 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2psidxideta_map, and libMesh::libmesh_assert().
|
inline |
Second derivatives of "xi" reference coordinate wrt physical coordinates.
Definition at line 381 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2xidxyz2_map, and libMesh::libmesh_assert().
Referenced by libMesh::HCurlFETransformation< OutputShape >::init_map_d2phi(), libMesh::HDivFETransformation< OutputShape >::init_map_d2phi(), libMesh::H1FETransformation< OutputShape >::init_map_d2phi(), and libMesh::H1FETransformation< OutputShape >::map_d2phi().
|
inline |
Definition at line 271 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2xyzdeta2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 299 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2xyzdetadzeta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 264 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2xyzdxi2_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 285 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2xyzdxideta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 292 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2xyzdxidzeta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 278 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2xyzdzeta2_map, and libMesh::libmesh_assert().
|
inline |
Second derivatives of "zeta" reference coordinate wrt physical coordinates.
Definition at line 395 of file fe_map.h.
References calculate_d2xyz, calculations_started, d2zetadxyz2_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_d2phi().
|
inline |
Definition at line 333 of file fe_map.h.
References calculate_dxyz, calculations_started, detadx_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 341 of file fe_map.h.
References calculate_dxyz, calculations_started, detady_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 349 of file fe_map.h.
References calculate_dxyz, calculations_started, detadz_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 423 of file fe_map.h.
References calculate_dxyz, calculations_started, dphideta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 566 of file fe_map.h.
References calculate_dxyz, calculations_started, dphideta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 416 of file fe_map.h.
References calculate_dxyz, calculations_started, dphidxi_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 559 of file fe_map.h.
References calculate_dxyz, calculations_started, dphidxi_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 430 of file fe_map.h.
References calculate_dxyz, calculations_started, dphidzeta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 573 of file fe_map.h.
References calculate_dxyz, calculations_started, dphidzeta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 494 of file fe_map.h.
References calculate_dxyz, calculations_started, dpsideta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 498 of file fe_map.h.
References calculate_dxyz, calculations_started, dpsideta_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 483 of file fe_map.h.
References calculate_dxyz, calculations_started, dpsidxi_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 487 of file fe_map.h.
References calculate_dxyz, calculations_started, dpsidxi_map, and libMesh::libmesh_assert().
|
inline |
Definition at line 309 of file fe_map.h.
References calculate_dxyz, calculations_started, dxidx_map, and libMesh::libmesh_assert().
Referenced by libMesh::HCurlFETransformation< OutputShape >::init_map_d2phi(), libMesh::HDivFETransformation< OutputShape >::init_map_d2phi(), libMesh::H1FETransformation< OutputShape >::init_map_d2phi(), libMesh::HCurlFETransformation< OutputShape >::init_map_dphi(), libMesh::HDivFETransformation< OutputShape >::init_map_dphi(), libMesh::H1FETransformation< OutputShape >::init_map_dphi(), libMesh::HDivFETransformation< OutputShape >::init_map_phi(), libMesh::HCurlFETransformation< OutputShape >::init_map_phi(), libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 317 of file fe_map.h.
References calculate_dxyz, calculations_started, dxidy_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 325 of file fe_map.h.
References calculate_dxyz, calculations_started, dxidz_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 247 of file fe_map.h.
References calculate_dxyz, calculations_started, dxyzdeta_map, and libMesh::libmesh_assert().
Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl(), and libMesh::HDivFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 239 of file fe_map.h.
References calculate_dxyz, calculations_started, dxyzdxi_map, and libMesh::libmesh_assert().
Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl(), and libMesh::HDivFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 255 of file fe_map.h.
References calculate_dxyz, calculations_started, dxyzdzeta_map, and libMesh::libmesh_assert().
Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl(), and libMesh::HDivFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 357 of file fe_map.h.
References calculate_dxyz, calculations_started, dzetadx_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 365 of file fe_map.h.
References calculate_dxyz, calculations_started, dzetady_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 373 of file fe_map.h.
References calculate_dxyz, calculations_started, dzetadz_map, and libMesh::libmesh_assert().
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 223 of file fe_map.h.
References calculate_dxyz, calculations_started, jac, and libMesh::libmesh_assert().
Referenced by libMesh::Elem::has_invertible_map(), libMesh::HCurlFETransformation< OutputShape >::map_curl(), libMesh::HDivFETransformation< OutputShape >::map_div(), and libMesh::HDivFETransformation< OutputShape >::map_phi().
|
inline |
Definition at line 231 of file fe_map.h.
References calculate_dxyz, calculations_started, JxW, and libMesh::libmesh_assert().
|
inline |
Definition at line 627 of file fe_map.h.
References calculate_dxyz, calculations_started, JxW, and libMesh::libmesh_assert().
|
inline |
Definition at line 444 of file fe_map.h.
References calculate_dxyz, calculations_started, libMesh::libmesh_assert(), and normals.
|
inline |
Definition at line 409 of file fe_map.h.
References calculate_xyz, calculations_started, libMesh::libmesh_assert(), and phi_map.
|
inline |
Definition at line 552 of file fe_map.h.
References calculate_xyz, calculations_started, libMesh::libmesh_assert(), and phi_map.
|
inline |
Definition at line 403 of file fe_map.h.
References psi_map.
|
inline |
Definition at line 477 of file fe_map.h.
References psi_map.
|
inline |
Definition at line 437 of file fe_map.h.
References calculate_dxyz, calculations_started, libMesh::libmesh_assert(), and tangents.
|
inline |
xyz
spatial locations of the quadrature points on the element. Definition at line 216 of file fe_map.h.
References calculate_xyz, calculations_started, libMesh::libmesh_assert(), and xyz.
void libMesh::FEMap::init_edge_shape_functions | ( | const std::vector< Point > & | qp, |
const Elem * | edge | ||
) |
Same as before, but for an edge.
This is used for some projection operators.
Definition at line 565 of file fe_boundary.C.
References calculate_d2xyz, calculate_dxyz, calculate_xyz, d2psidxi2_map, libMesh::Elem::default_order(), determine_calculations(), dpsidxi_map, libMesh::libmesh_assert(), map_fe_type(), libMesh::FEInterface::n_shape_functions(), psi_map, libMesh::FEInterface::shape_deriv_function(), libMesh::FEInterface::shape_function(), and libMesh::FEInterface::shape_second_deriv_function().
void libMesh::FEMap::init_face_shape_functions | ( | const std::vector< Point > & | qp, |
const Elem * | side | ||
) |
Initializes the reference to physical element map for a side.
This is used for boundary integration.
Definition at line 437 of file fe_boundary.C.
References calculate_d2xyz, calculate_dxyz, calculate_xyz, d2psideta2_map, d2psidxi2_map, d2psidxideta_map, libMesh::Elem::default_order(), determine_calculations(), dpsideta_map, dpsidxi_map, libMesh::libmesh_assert(), map_fe_type(), libMesh::FEInterface::n_shape_functions(), psi_map, libMesh::FEInterface::shape_deriv_function(), libMesh::FEInterface::shape_function(), and libMesh::FEInterface::shape_second_deriv_function().
void libMesh::FEMap::init_reference_to_physical_map | ( | const std::vector< Point > & | qp, |
const Elem * | elem | ||
) |
Definition at line 108 of file fe_map.C.
References libMesh::FEInterface::all_shape_derivs(), libMesh::FEInterface::all_shapes(), calculate_d2xyz, calculate_dxyz, calculate_xyz, d2phideta2_map, d2phidetadzeta_map, d2phidxi2_map, d2phidxideta_map, d2phidxidzeta_map, d2phidzeta2_map, libMesh::Elem::default_order(), determine_calculations(), dphideta_map, dphidxi_map, dphidzeta_map, libMesh::Elem::is_linear(), map_fe_type(), libMesh::FEInterface::n_shape_functions(), phi_map, libMesh::FEInterface::shape_deriv_function(), and libMesh::FEInterface::shape_second_deriv_function().
Referenced by libMesh::Elem::has_invertible_map().
|
inline |
|
static |
p
located in physical space. This function requires inverting the (possibly 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}} \)When secure == true, the following checks are enabled:
In DEBUG mode only: .) dim==1,2: throw an error if det(J) <= 0 for any Newton iteration. .) Print warning for every iteration beyond max_cnt in which the Newton scheme has not converged.
In !DEBUG mode only: .) Print a single warning (1 warning for the entire simulation) if the Newton scheme ever requires more than max_cnt iterations.
In both DEBUG and !DEBUG modes: .) dim==3: Throw an exception for singular Jacobian. .) Throw an error if the Newton iteration has not converged in 2*max_cnt iterations.
In addition to the checks above, the "extra_checks" parameter can be used to turn on some additional tests. In particular, when extra_checks == true and compiled in DEBUG mode: .) Print a warning if p != map(inverse_map(p)) to within tolerance. .) Print a warning if the inverse-mapped point is not on the reference element to within tolerance.
Definition at line 1626 of file fe_map.C.
References libMesh::TypeVector< T >::add(), dim, libMesh::err, libMesh::DofObject::id(), libMesh::Elem::infinite(), libMesh::invalid_uint, libMesh::InfFEMap::inverse_map(), libMesh::libmesh_assert(), libMesh::libmesh_ignore(), libMesh::Elem::local_singular_node(), map(), map_deriv(), libMesh::Elem::master_point(), libMesh::Elem::n_nodes(), libMesh::TypeVector< T >::norm(), libMesh::FEAbstract::on_reference_element(), libMesh::Elem::print_info(), libMesh::Real, and libMesh::Elem::type().
Referenced by libMesh::MeshFunction::_gradient_on_elem(), libMesh::HPCoarsenTest::add_projection(), assemble_ellipticdg(), assemble_poisson(), InfFERadialTest::base_point(), libMesh::FEMContext::build_new_fe(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), compute_face_map(), compute_jacobian(), 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(), compute_residual(), libMesh::MeshFunction::discontinuous_value(), libMesh::DTKEvaluator::evaluate(), libMesh::MeshFunction::hessian(), libMesh::InfFE< Dim, T_radial, T_map >::inf_compute_constraints(), libMesh::RBEIMConstruction::initialize_qp_data(), libMesh::InfFEMap::inverse_map(), inverse_map(), libMesh::FE< Dim, LAGRANGE_VEC >::inverse_map(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), libMesh::MeshFunction::operator()(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectEdges::operator()(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectSides::operator()(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectInteriors::operator()(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::Elem::point_test(), libMesh::System::point_value(), libMesh::JumpErrorEstimator::reinit_sides(), libMesh::HPCoarsenTest::select_refinement(), NavierSystem::side_constraint(), FETest< order, family, elem_type >::testLoop(), and MeshInputTest::testMasterCenters().
|
static |
Takes a number points in physical space (in the physical_points
vector) and finds their location on the reference element for the input element elem
.
The values on the reference element are returned in the vector reference_points
. The other parameters have the same meaning as the single Point version of inverse_map() above.
Definition at line 2032 of file fe_map.C.
References dim, libMesh::Elem::infinite(), libMesh::InfFEMap::inverse_map(), inverse_map(), and libMesh::libmesh_ignore().
|
static |
p
located on the reference element. Definition at line 2067 of file fe_map.C.
References libMesh::TypeVector< T >::add_scaled(), libMesh::Elem::default_order(), dim, libMesh::Elem::infinite(), libMesh::libmesh_assert(), libMesh::libmesh_ignore(), libMesh::InfFEMap::map(), map_fe_type(), libMesh::FEInterface::n_shape_functions(), libMesh::Elem::point(), and libMesh::FEInterface::shape_function().
Referenced by libMesh::RBEIMConstruction::initialize_qp_data(), inverse_map(), libMesh::InfFEMap::map(), libMesh::FE< Dim, LAGRANGE_VEC >::map(), libMesh::Elem::point_test(), and MeshInputTest::testMasterCenters().
|
static |
j
of d(xyz)/d(xi eta zeta) (in physical space) of the point p
located on the reference element. Definition at line 2101 of file fe_map.C.
References libMesh::TypeVector< T >::add_scaled(), libMesh::Elem::default_order(), dim, libMesh::Elem::infinite(), libMesh::libmesh_assert(), libMesh::libmesh_ignore(), map_fe_type(), libMesh::FEInterface::n_shape_functions(), libMesh::Elem::point(), and libMesh::FEInterface::shape_deriv_function().
Referenced by compute_face_map(), inverse_map(), libMesh::FE< Dim, LAGRANGE_VEC >::map_eta(), libMesh::FE< Dim, LAGRANGE_VEC >::map_xi(), and libMesh::FE< Dim, LAGRANGE_VEC >::map_zeta().
Definition at line 45 of file fe_map.C.
References libMesh::LAGRANGE, libMesh::LAGRANGE_MAP, libMesh::Elem::mapping_type(), libMesh::RATIONAL_BERNSTEIN, and libMesh::RATIONAL_BERNSTEIN_MAP.
Referenced by libMesh::FEMContext::_do_elem_position_set(), compute_face_map(), libMesh::FEAbstract::compute_node_constraints(), libMesh::FEAbstract::compute_periodic_node_constraints(), libMesh::FEMContext::elem_position_get(), init_edge_shape_functions(), init_face_shape_functions(), init_reference_to_physical_map(), map(), map_deriv(), libMesh::Elem::true_centroid(), and libMesh::Elem::volume().
void libMesh::FEMap::print_JxW | ( | std::ostream & | os | ) | const |
Prints the Jacobian times the weight for each quadrature point.
Definition at line 1490 of file fe_map.C.
References libMesh::index_range(), and JxW.
void libMesh::FEMap::print_xyz | ( | std::ostream & | os | ) | const |
Prints the spatial location of each quadrature point (on the physical element).
Definition at line 1498 of file fe_map.C.
References libMesh::index_range(), and xyz.
|
protected |
A utility function for use by compute_*_map.
Definition at line 1200 of file fe_map.C.
References calculate_d2xyz, calculate_dxyz, calculate_xyz, d2etadxyz2_map, d2xidxyz2_map, d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, d2zetadxyz2_map, detadx_map, detady_map, detadz_map, determine_calculations(), dim, dxidx_map, dxidy_map, dxidz_map, dxyzdeta_map, dxyzdxi_map, dxyzdzeta_map, dzetadx_map, dzetady_map, dzetadz_map, libMesh::index_range(), jac, JxW, and xyz.
Referenced by compute_affine_map(), compute_map(), and compute_null_map().
|
inline |
Set the Jacobian tolerance used for determining when the mapping fails.
The mapping is determined to fail if jac <= jacobian_tolerance.
Definition at line 641 of file fe_map.h.
References jacobian_tolerance.
|
private |
Work vector for compute_affine_map()
Definition at line 1044 of file fe_map.h.
Referenced by compute_affine_map(), and compute_map().
|
mutableprotected |
Should we calculate mapping hessians?
Definition at line 1023 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), determine_calculations(), get_curvatures(), get_d2etadxyz2(), get_d2phideta2_map(), get_d2phidetadzeta_map(), get_d2phidxi2_map(), get_d2phidxideta_map(), get_d2phidxidzeta_map(), get_d2phidzeta2_map(), get_d2psideta2(), get_d2psidxi2(), get_d2psidxideta(), get_d2xidxyz2(), get_d2xyzdeta2(), get_d2xyzdetadzeta(), get_d2xyzdxi2(), get_d2xyzdxideta(), get_d2xyzdxidzeta(), get_d2xyzdzeta2(), get_d2zetadxyz2(), init_edge_shape_functions(), init_face_shape_functions(), init_reference_to_physical_map(), and resize_quadrature_map_vectors().
|
mutableprotected |
Should we calculate mapping gradients?
Definition at line 1016 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), determine_calculations(), get_detadx(), get_detady(), get_detadz(), get_dphideta_map(), get_dphidxi_map(), get_dphidzeta_map(), get_dpsideta(), get_dpsidxi(), get_dxidx(), get_dxidy(), get_dxidz(), get_dxyzdeta(), get_dxyzdxi(), get_dxyzdzeta(), get_dzetadx(), get_dzetady(), get_dzetadz(), get_jacobian(), get_JxW(), get_normals(), get_tangents(), init_edge_shape_functions(), init_face_shape_functions(), init_reference_to_physical_map(), and resize_quadrature_map_vectors().
|
mutableprotected |
Should we calculate physical point locations?
Definition at line 1011 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), libMesh::FEXYZMap::FEXYZMap(), get_phi_map(), get_xyz(), init_edge_shape_functions(), init_face_shape_functions(), init_reference_to_physical_map(), and resize_quadrature_map_vectors().
|
mutableprotected |
Have calculations with this object already been started? Then all get_* functions should already have been called.
Definition at line 1006 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), determine_calculations(), get_curvatures(), get_d2etadxyz2(), get_d2phideta2_map(), get_d2phidetadzeta_map(), get_d2phidxi2_map(), get_d2phidxideta_map(), get_d2phidxidzeta_map(), get_d2phidzeta2_map(), get_d2psideta2(), get_d2psidxi2(), get_d2psidxideta(), get_d2xidxyz2(), get_d2xyzdeta2(), get_d2xyzdetadzeta(), get_d2xyzdxi2(), get_d2xyzdxideta(), get_d2xyzdxidzeta(), get_d2xyzdzeta2(), get_d2zetadxyz2(), get_detadx(), get_detady(), get_detadz(), get_dphideta_map(), get_dphidxi_map(), get_dphidzeta_map(), get_dpsideta(), get_dpsidxi(), get_dxidx(), get_dxidy(), get_dxidz(), get_dxyzdeta(), get_dxyzdxi(), get_dxyzdzeta(), get_dzetadx(), get_dzetady(), get_dzetadz(), get_jacobian(), get_JxW(), get_normals(), get_phi_map(), get_tangents(), and get_xyz().
|
protected |
The mean curvature (= one half the sum of the principal curvatures) on the boundary at the quadrature points.
The mean curvature is a scalar value.
Definition at line 988 of file fe_map.h.
Referenced by compute_edge_map(), compute_face_map(), and get_curvatures().
|
protected |
Second derivatives of "eta" reference coordinate wrt physical coordinates.
At each qp: (eta_{xx}, eta_{xy}, eta_{xz}, eta_{yy}, eta_{yz}, eta_{zz})
Definition at line 867 of file fe_map.h.
Referenced by compute_inverse_map_second_derivs(), compute_single_point_map(), get_d2etadxyz2(), and resize_quadrature_map_vectors().
|
protected |
Map for the second derivative, d^2(phi)/d(eta)^2.
Definition at line 916 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_d2phideta2_map(), and init_reference_to_physical_map().
|
protected |
Map for the second derivative, d^2(phi)/d(eta)d(zeta).
Definition at line 921 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_d2phidetadzeta_map(), and init_reference_to_physical_map().
|
protected |
Map for the second derivative, d^2(phi)/d(xi)^2.
Definition at line 901 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_d2phidxi2_map(), and init_reference_to_physical_map().
|
protected |
Map for the second derivative, d^2(phi)/d(xi)d(eta).
Definition at line 906 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_d2phidxideta_map(), and init_reference_to_physical_map().
|
protected |
Map for the second derivative, d^2(phi)/d(xi)d(zeta).
Definition at line 911 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_d2phidxidzeta_map(), and init_reference_to_physical_map().
|
protected |
Map for the second derivative, d^2(phi)/d(zeta)^2.
Definition at line 926 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_d2phidzeta2_map(), and init_reference_to_physical_map().
|
protected |
Map for the second derivatives (in eta) of the side shape functions.
Useful for computing the curvature at the quadrature points.
Definition at line 968 of file fe_map.h.
Referenced by compute_face_map(), get_d2psideta2(), and init_face_shape_functions().
|
protected |
Map for the second derivatives (in xi) of the side shape functions.
Useful for computing the curvature at the quadrature points.
Definition at line 954 of file fe_map.h.
Referenced by compute_edge_map(), compute_face_map(), get_d2psidxi2(), init_edge_shape_functions(), and init_face_shape_functions().
|
protected |
Map for the second (cross) derivatives in xi, eta of the side shape functions.
Useful for computing the curvature at the quadrature points.
Definition at line 961 of file fe_map.h.
Referenced by compute_face_map(), get_d2psidxideta(), and init_face_shape_functions().
|
protected |
Second derivatives of "xi" reference coordinate wrt physical coordinates.
At each qp: (xi_{xx}, xi_{xy}, xi_{xz}, xi_{yy}, xi_{yz}, xi_{zz})
Definition at line 861 of file fe_map.h.
Referenced by compute_inverse_map_second_derivs(), compute_single_point_map(), get_d2xidxyz2(), and resize_quadrature_map_vectors().
|
protected |
Vector of second partial derivatives in eta: d^2(x)/d(eta)^2.
Definition at line 778 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_d2xyzdeta2(), and resize_quadrature_map_vectors().
|
protected |
Vector of mixed second partial derivatives in eta-zeta: d^2(x)/d(eta)d(zeta) d^2(y)/d(eta)d(zeta) d^2(z)/d(eta)d(zeta)
Definition at line 790 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_d2xyzdetadzeta(), and resize_quadrature_map_vectors().
|
protected |
Vector of second partial derivatives in xi: d^2(x)/d(xi)^2, d^2(y)/d(xi)^2, d^2(z)/d(xi)^2.
Definition at line 766 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_d2xyzdxi2(), and resize_quadrature_map_vectors().
|
protected |
Vector of mixed second partial derivatives in xi-eta: d^2(x)/d(xi)d(eta) d^2(y)/d(xi)d(eta) d^2(z)/d(xi)d(eta)
Definition at line 772 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_d2xyzdxideta(), and resize_quadrature_map_vectors().
|
protected |
Vector of second partial derivatives in xi-zeta: d^2(x)/d(xi)d(zeta), d^2(y)/d(xi)d(zeta), d^2(z)/d(xi)d(zeta)
Definition at line 784 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_d2xyzdxidzeta(), and resize_quadrature_map_vectors().
|
protected |
Vector of second partial derivatives in zeta: d^2(x)/d(zeta)^2.
Definition at line 796 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_d2xyzdzeta2(), and resize_quadrature_map_vectors().
|
protected |
Second derivatives of "zeta" reference coordinate wrt physical coordinates.
At each qp: (zeta_{xx}, zeta_{xy}, zeta_{xz}, zeta_{yy}, zeta_{yz}, zeta_{zz})
Definition at line 873 of file fe_map.h.
Referenced by compute_inverse_map_second_derivs(), compute_single_point_map(), get_d2zetadxyz2(), and resize_quadrature_map_vectors().
|
protected |
Map for partial derivatives: d(eta)/d(x).
Needed for the Jacobian.
Definition at line 823 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_detadx(), and resize_quadrature_map_vectors().
|
protected |
Map for partial derivatives: d(eta)/d(y).
Needed for the Jacobian.
Definition at line 829 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_detady(), and resize_quadrature_map_vectors().
|
protected |
Map for partial derivatives: d(eta)/d(z).
Needed for the Jacobian.
Definition at line 835 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_detadz(), and resize_quadrature_map_vectors().
|
protected |
Map for the derivative, d(phi)/d(eta).
Definition at line 889 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_dphideta_map(), and init_reference_to_physical_map().
|
protected |
Map for the derivative, d(phi)/d(xi).
Definition at line 884 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_dphidxi_map(), and init_reference_to_physical_map().
|
protected |
Map for the derivative, d(phi)/d(zeta).
Definition at line 894 of file fe_map.h.
Referenced by add_calculations(), compute_single_point_map(), get_dphidzeta_map(), and init_reference_to_physical_map().
|
protected |
Map for the derivative of the side function, d(psi)/d(eta).
Definition at line 945 of file fe_map.h.
Referenced by compute_face_map(), get_dpsideta(), and init_face_shape_functions().
|
protected |
Map for the derivative of the side functions, d(psi)/d(xi).
Definition at line 939 of file fe_map.h.
Referenced by compute_edge_map(), compute_face_map(), get_dpsidxi(), init_edge_shape_functions(), and init_face_shape_functions().
|
protected |
Map for partial derivatives: d(xi)/d(x).
Needed for the Jacobian.
Definition at line 804 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_dxidx(), and resize_quadrature_map_vectors().
|
protected |
Map for partial derivatives: d(xi)/d(y).
Needed for the Jacobian.
Definition at line 810 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_dxidy(), and resize_quadrature_map_vectors().
|
protected |
Map for partial derivatives: d(xi)/d(z).
Needed for the Jacobian.
Definition at line 816 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_dxidz(), and resize_quadrature_map_vectors().
|
protected |
Vector of partial derivatives: d(x)/d(eta), d(y)/d(eta), d(z)/d(eta)
Definition at line 752 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), dxdeta_map(), dydeta_map(), dzdeta_map(), get_dxyzdeta(), and resize_quadrature_map_vectors().
|
protected |
Vector of partial derivatives: d(x)/d(xi), d(y)/d(xi), d(z)/d(xi)
Definition at line 746 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), dxdxi_map(), dydxi_map(), dzdxi_map(), get_dxyzdxi(), and resize_quadrature_map_vectors().
|
protected |
Vector of partial derivatives: d(x)/d(zeta), d(y)/d(zeta), d(z)/d(zeta)
Definition at line 758 of file fe_map.h.
Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), dxdzeta_map(), dydzeta_map(), dzdzeta_map(), get_dxyzdzeta(), and resize_quadrature_map_vectors().
|
protected |
Map for partial derivatives: d(zeta)/d(x).
Needed for the Jacobian.
Definition at line 842 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_dzetadx(), and resize_quadrature_map_vectors().
|
protected |
Map for partial derivatives: d(zeta)/d(y).
Needed for the Jacobian.
Definition at line 848 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_dzetady(), and resize_quadrature_map_vectors().
|
protected |
Map for partial derivatives: d(zeta)/d(z).
Needed for the Jacobian.
Definition at line 854 of file fe_map.h.
Referenced by compute_affine_map(), compute_inverse_map_second_derivs(), compute_null_map(), compute_single_point_map(), get_dzetadz(), and resize_quadrature_map_vectors().
|
protected |
Jacobian values at quadrature points.
Definition at line 995 of file fe_map.h.
Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_jacobian(), and resize_quadrature_map_vectors().
|
protected |
The Jacobian tolerance used for determining when the mapping fails.
The mapping is determined to fail if jac <= jacobian_tolerance. If not set by the user, this number defaults to 0
Definition at line 1032 of file fe_map.h.
Referenced by compute_single_point_map(), and set_jacobian_tolerance().
|
protected |
Jacobian*Weight values at quadrature points.
Definition at line 1000 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), get_JxW(), print_JxW(), and resize_quadrature_map_vectors().
|
protected |
Normal vectors on boundary at quadrature points.
Definition at line 980 of file fe_map.h.
Referenced by compute_edge_map(), compute_face_map(), and get_normals().
|
protected |
Map for the shape function phi.
Definition at line 879 of file fe_map.h.
Referenced by add_calculations(), compute_affine_map(), compute_single_point_map(), get_phi_map(), and init_reference_to_physical_map().
|
protected |
Map for the side shape functions, psi.
Definition at line 933 of file fe_map.h.
Referenced by compute_edge_map(), compute_face_map(), get_psi(), init_edge_shape_functions(), and init_face_shape_functions().
|
protected |
Tangent vectors on boundary at quadrature points.
Definition at line 975 of file fe_map.h.
Referenced by compute_edge_map(), compute_face_map(), and get_tangents().
|
protected |
The spatial locations of the quadrature points.
Definition at line 740 of file fe_map.h.
Referenced by compute_affine_map(), compute_edge_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), get_xyz(), print_xyz(), and resize_quadrature_map_vectors().