Go to the documentation of this file.
   18 #include "libmesh/fe_xyz_map.h" 
   19 #include "libmesh/elem.h" 
   20 #include "libmesh/libmesh_logging.h" 
   29   LOG_SCOPE(
"compute_face_map()", 
"FEXYZMap");
 
   32   const unsigned int n_qp = cast_int<unsigned int>(qw.size());
 
   41           this->
xyz.resize(n_qp);
 
   43 #ifdef  LIBMESH_ENABLE_SECOND_DERIVATIVES 
   50           this->
JxW.resize(n_qp);
 
   55         for (
unsigned int p=0; p<n_qp; p++)
 
   57             this->
tangents[p].resize(LIBMESH_DIM-1); 
 
   60 #ifdef  LIBMESH_ENABLE_SECOND_DERIVATIVES 
   66         for (
unsigned int i=0,
 
   67              n_psi_map = cast_int<unsigned int>(this->
psi_map.size());
 
   72             for (
unsigned int p=0; p<n_qp; p++) 
 
   74                 this->
xyz[p].add_scaled          (side_point, this->
psi_map[i][p]);
 
   76 #ifdef  LIBMESH_ENABLE_SECOND_DERIVATIVES 
   83         for (
unsigned int p=0; p<n_qp; p++)
 
   89 #if LIBMESH_DIM == 3  // Only good in 3D space 
   92 #ifdef  LIBMESH_ENABLE_SECOND_DERIVATIVES 
  104             libmesh_assert_not_equal_to (denominator, 0);
 
  105             this->
curvatures[p] = numerator / denominator;
 
  111         for (
unsigned int p=0; p<n_qp; p++)
 
  116             libmesh_assert_greater (the_jac, 0.);
 
  118             this->
JxW[p] = the_jac*qw[p];
 
  128           this->
xyz.resize(n_qp);
 
  131 #ifdef  LIBMESH_ENABLE_SECOND_DERIVATIVES 
  140           this->
JxW.resize(n_qp);
 
  144         for (
unsigned int p=0; p<n_qp; p++)
 
  146             this->
tangents[p].resize(LIBMESH_DIM-1); 
 
  150 #ifdef  LIBMESH_ENABLE_SECOND_DERIVATIVES 
  158         for (
unsigned int i=0,
 
  159              n_psi_map = cast_int<unsigned int>(this->
psi_map.size());
 
  164             for (
unsigned int p=0; p<n_qp; p++) 
 
  166                 this->
xyz[p].add_scaled             (side_point, this->
psi_map[i][p]);
 
  169 #ifdef  LIBMESH_ENABLE_SECOND_DERIVATIVES 
  178         for (
unsigned int p=0; p<n_qp; p++)
 
  185 #ifdef  LIBMESH_ENABLE_SECOND_DERIVATIVES 
  197             const Real G  =  this->dxyzdeta_map[p].norm_sq();
 
  199             const Real numerator   = E*N -2.*F*M + G*L;
 
  200             const Real denominator = E*G - F*F;
 
  201             libmesh_assert_not_equal_to (denominator, 0.);
 
  202             this->
curvatures[p] = 0.5*numerator/denominator;
 
  208         for (
unsigned int p=0; p<n_qp; p++)
 
  218             const Real g21 = g12;
 
  227             libmesh_assert_greater (the_jac, 0.);
 
  229             this->
JxW[p] = the_jac*qw[p];
 
  235       libmesh_error_msg(
"Invalid dim = " << 
dim);
 
  
std::vector< std::vector< Real > > d2psidxi2_map
Map for the second derivatives (in xi) of the side shape functions.
 
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 > > dpsidxi_map
Map for the derivative of the side functions, d(psi)/d(xi).
 
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.
 
virtual void compute_face_map(int dim, const std::vector< Real > &qw, const Elem *side) override
Special implementation for XYZ finite elements.
 
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.
 
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
 
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< Point > normals
Normal vectors on boundary at quadrature points.
 
std::vector< Point > xyz
The spatial locations of the quadrature points.
 
const Point & point(const unsigned int i) const
 
std::vector< std::vector< Real > > dpsideta_map
Map for the derivative of the side function, d(psi)/d(eta).
 
std::vector< std::vector< Point > > tangents
Tangent vectors on boundary at quadrature points.
 
std::vector< std::vector< Real > > psi_map
Map for the side shape functions, psi.
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
std::vector< std::vector< Real > > d2psideta2_map
Map for the second derivatives (in eta) of the side shape functions.
 
Real dzdeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
 
TypeVector< T > unit() const
 
This is the base class from which all geometric element types are derived.
 
Real dzdxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
 
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
 
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)
 
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(...
 
std::vector< RealGradient > d2xyzdeta2_map
Vector of second partial derivatives in eta: d^2(x)/d(eta)^2.
 
std::vector< std::vector< Real > > d2psidxideta_map
Map for the second (cross) derivatives in xi, eta of the side shape functions.
 
std::vector< Real > curvatures
The mean curvature (= one half the sum of the principal curvatures) on the boundary at the quadrature...