Go to the documentation of this file.
20 #ifndef LIBMESH_FE_MAP_H
21 #define LIBMESH_FE_MAP_H
24 #include "libmesh/reference_counted_object.h"
25 #include "libmesh/point.h"
26 #include "libmesh/vector_value.h"
27 #include "libmesh/fe_type.h"
54 static std::unique_ptr<FEMap>
build(
FEType fe_type);
58 template<
unsigned int Dim>
70 const std::vector<Real> & qw,
73 const std::vector<const Node *> & elem_nodes,
74 bool compute_second_derivatives);
83 const std::vector<Real> & qw,
92 const std::vector<Real> & qw);
104 const std::vector<Real> & qw,
106 bool calculate_d2phi);
112 const std::vector<Real> & qw,
119 const std::vector<Real> & qw,
126 template<
unsigned int Dim>
134 template<
unsigned int Dim>
144 const Point & reference_point);
152 const unsigned int j,
153 const Point & reference_point);
170 const bool secure =
true);
186 const std::vector<Point> & physical_points,
187 std::vector<Point> & reference_points,
189 const bool secure =
true);
238 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
356 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
382 const std::vector<std::vector<Real>> &
get_psi()
const
427 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
481 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
526 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
556 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
624 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
732 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
770 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
828 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
868 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
919 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
954 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
990 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1003 template <
unsigned int Dim, FEFamily T>
1028 #endif // LIBMESH_FE_MAP_H
const std::vector< Real > & get_jacobian() const
std::vector< std::vector< Real > > d2psidxi2_map
Map for the second derivatives (in xi) of the side shape functions.
const std::vector< std::vector< Real > > & get_psi() const
std::vector< Real > dzetadx_map
Map for partial derivatives: d(zeta)/d(x).
bool calculate_d2xyz
Should we calculate mapping hessians?
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...
Real dydeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
std::vector< std::vector< Real > > & get_d2psidxideta()
const std::vector< std::vector< Real > > & get_d2zetadxyz2() const
Second derivatives of "zeta" reference coordinate wrt physical coordinates.
const std::vector< RealGradient > & get_d2xyzdzeta2() const
std::vector< std::vector< Real > > & get_d2phidxidzeta_map()
const std::vector< RealGradient > & get_d2xyzdxi2() const
void set_jacobian_tolerance(Real tol)
Set the Jacobian tolerance used for determining when the mapping fails.
const std::vector< RealGradient > & get_d2xyzdeta2() const
const std::vector< Real > & get_dzetadx() const
const std::vector< Real > & get_dxidx() const
std::vector< std::vector< Real > > dpsidxi_map
Map for the derivative of the side functions, d(psi)/d(xi).
std::vector< std::vector< Real > > & get_d2phidzeta2_map()
static Point map_deriv(const unsigned int dim, const Elem *elem, const unsigned int j, const Point &reference_point)
std::vector< Real > dzetadz_map
Map for partial derivatives: d(zeta)/d(z).
const std::vector< RealGradient > & get_dxyzdxi() const
std::vector< Real > dxidy_map
Map for partial derivatives: d(xi)/d(y).
const std::vector< std::vector< Real > > & get_d2psidxideta() const
const std::vector< std::vector< Real > > & get_dphidxi_map() const
bool calculate_xyz
Should we calculate physical point locations?
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.
The libMesh namespace provides an interface to certain functionality in the library.
const std::vector< Real > & get_dzetadz() const
std::vector< Real > jac
Jacobian values at quadrature points.
const std::vector< Real > & get_dzetady() const
Real dxdzeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
std::vector< std::vector< Real > > d2phidxideta_map
Map for the second derivative, d^2(phi)/d(xi)d(eta).
std::vector< Real > JxW
Jacobian*Weight values at quadrature points.
Real dxdxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
static const Real TOLERANCE
std::vector< std::vector< Real > > & get_psi()
void print_JxW(std::ostream &os) const
Prints the Jacobian times the weight for each quadrature point.
static Point map(const unsigned int dim, const Elem *elem, const Point &reference_point)
const std::vector< Real > & get_JxW() const
std::vector< const Node * > _elem_nodes
Work vector for compute_affine_map()
Real dydzeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
const std::vector< Point > & get_xyz() const
std::vector< RealGradient > dxyzdeta_map
Vector of partial derivatives: d(x)/d(eta), d(y)/d(eta), d(z)/d(eta)
Real dydxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
std::vector< std::vector< Real > > & get_d2psidxi2()
const std::vector< std::vector< Real > > & get_d2xidxyz2() const
Second derivatives of "xi" reference coordinate wrt physical coordinates.
std::vector< Point > normals
Normal vectors on boundary at quadrature points.
const std::vector< Real > & get_dxidy() const
std::vector< Point > xyz
The spatial locations of the quadrature points.
void init_face_shape_functions(const std::vector< Point > &qp, const Elem *side)
Initializes the reference to physical element map for a side.
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.
std::vector< std::vector< Real > > & get_dphideta_map()
std::vector< std::vector< Real > > d2phidxi2_map
Map for the second derivative, d^2(phi)/d(xi)^2.
virtual void compute_null_map(const unsigned int dim, const std::vector< Real > &qw)
Assign a fake jacobian and some other additional data fields.
std::vector< Real > detady_map
Map for partial derivatives: d(eta)/d(y).
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.
void print_xyz(std::ostream &os) const
Prints the spatial location of each quadrature point (on the physical element).
const std::vector< RealGradient > & get_dxyzdeta() const
virtual void compute_face_map(int dim, const std::vector< Real > &qw, const Elem *side)
Same as compute_map, but for a side.
const std::vector< std::vector< Real > > & get_d2psidxi2() const
const std::vector< Real > & get_detady() const
const std::vector< RealGradient > & get_d2xyzdxidzeta() const
void determine_calculations()
Determine which values are to be calculated.
std::vector< std::vector< Real > > & get_phi_map()
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),...
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
void resize_quadrature_map_vectors(const unsigned int dim, unsigned int n_qp)
A utility function for use by compute_*_map.
std::vector< std::vector< Real > > dpsideta_map
Map for the derivative of the side function, d(psi)/d(eta).
std::vector< std::vector< Real > > & get_dpsideta()
std::vector< Real > detadz_map
Map for partial derivatives: d(eta)/d(z).
std::vector< std::vector< Point > > tangents
Tangent vectors on boundary at quadrature points.
std::vector< std::vector< Real > > & get_d2phidxideta_map()
A specific instantiation of the FEBase class.
std::vector< std::vector< Real > > psi_map
Map for the side shape functions, psi.
std::vector< std::vector< Real > > d2phidzeta2_map
Map for the second derivative, d^2(phi)/d(zeta)^2.
A Point defines a location in LIBMESH_DIM dimensional Real space.
bool calculate_dxyz
Should we calculate mapping gradients?
const std::vector< RealGradient > & get_d2xyzdxideta() const
std::vector< RealGradient > d2xyzdzeta2_map
Vector of second partial derivatives in zeta: d^2(x)/d(zeta)^2.
std::vector< std::vector< Real > > dphideta_map
Map for the derivative, d(phi)/d(eta).
static std::unique_ptr< FEMap > build(FEType fe_type)
const std::vector< std::vector< Real > > & get_dphidzeta_map() const
const std::vector< Real > & get_detadx() const
std::vector< std::vector< Real > > d2xidxyz2_map
Second derivatives of "xi" reference coordinate wrt physical coordinates.
std::vector< std::vector< Real > > d2psideta2_map
Map for the second derivatives (in eta) of the side shape functions.
void compute_edge_map(int dim, const std::vector< Real > &qw, const Elem *side)
Same as before, but for an edge.
Real dzdeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
static FEFamily map_fe_type(const Elem &elem)
std::vector< std::vector< Real > > d2etadxyz2_map
Second derivatives of "eta" reference coordinate wrt physical coordinates.
std::vector< Real > detadx_map
Map for partial derivatives: d(eta)/d(x).
const std::vector< std::vector< Real > > & get_d2psideta2() const
void init_reference_to_physical_map(const std::vector< Point > &qp, const Elem *elem)
Real jacobian_tolerance
The Jacobian tolerance used for determining when the mapping fails.
void init_edge_shape_functions(const std::vector< Point > &qp, const Elem *edge)
Same as before, but for an edge.
std::vector< std::vector< Real > > dphidxi_map
Map for the derivative, d(phi)/d(xi).
std::vector< std::vector< Real > > & get_dphidzeta_map()
void compute_inverse_map_second_derivs(unsigned p)
A helper function used by FEMap::compute_single_point_map() to compute second derivatives of the inve...
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
std::vector< Real > & get_JxW()
std::vector< std::vector< Real > > d2phidetadzeta_map
Map for the second derivative, d^2(phi)/d(eta)d(zeta).
std::vector< std::vector< Real > > phi_map
Map for the shape function phi.
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.
std::vector< std::vector< Real > > & get_dpsidxi()
std::vector< std::vector< Real > > & get_d2phideta2_map()
std::vector< RealGradient > dxyzdzeta_map
Vector of partial derivatives: d(x)/d(zeta), d(y)/d(zeta), d(z)/d(zeta)
std::vector< Real > dzetady_map
Map for partial derivatives: d(zeta)/d(y).
std::vector< std::vector< Real > > dphidzeta_map
Map for the derivative, d(phi)/d(zeta).
std::vector< std::vector< Real > > & get_d2phidxi2_map()
This is the base class from which all geometric element types are derived.
const std::vector< Real > & get_curvatures() const
const std::vector< Point > & get_normals() const
const std::vector< Real > & get_detadz() const
const std::vector< std::vector< Real > > & get_phi_map() const
const std::vector< std::vector< Real > > & get_dpsidxi() const
std::vector< std::vector< Real > > d2phideta2_map
Map for the second derivative, d^2(phi)/d(eta)^2.
Class contained in FE that encapsulates mapping (i.e.
std::vector< std::vector< Real > > d2zetadxyz2_map
Second derivatives of "zeta" reference coordinate wrt physical coordinates.
std::vector< std::vector< Real > > & get_d2psideta2()
std::vector< std::vector< Real > > d2phidxidzeta_map
Map for the second derivative, d^2(phi)/d(xi)d(zeta).
const std::vector< std::vector< Real > > & get_dpsideta() const
const std::vector< RealGradient > & get_dxyzdzeta() const
Real dzdxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
Real dzdzeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real dxdeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
std::vector< RealGradient > dxyzdxi_map
Vector of partial derivatives: d(x)/d(xi), d(y)/d(xi), d(z)/d(xi)
const std::vector< std::vector< Point > > & get_tangents() const
std::vector< std::vector< Real > > & get_dphidxi_map()
std::vector< Real > dxidx_map
Map for partial derivatives: d(xi)/d(x).
const std::vector< Real > & get_dxidz() const
std::vector< std::vector< Real > > & get_d2phidetadzeta_map()
std::vector< Real > dxidz_map
Map for partial derivatives: d(xi)/d(z).
bool calculations_started
Have calculations with this object already been started? Then all get_* functions should already have...
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(...
const std::vector< std::vector< Real > > & get_dphideta_map() const
std::vector< RealGradient > d2xyzdeta2_map
Vector of second partial derivatives in eta: d^2(x)/d(eta)^2.
const std::vector< std::vector< Real > > & get_d2etadxyz2() const
Second derivatives of "eta" reference coordinate wrt physical coordinates.
std::vector< std::vector< Real > > d2psidxideta_map
Map for the second (cross) derivatives in xi, eta of the side shape functions.
const std::vector< RealGradient > & get_d2xyzdetadzeta() const
std::vector< Real > curvatures
The mean curvature (= one half the sum of the principal curvatures) on the boundary at the quadrature...