$extrastylesheet
libMesh::FEMap Class Reference

#include <fe_map.h>

Inheritance diagram for libMesh::FEMap:

List of all members.

Public Member Functions

 FEMap ()
virtual ~FEMap ()
template<unsigned int Dim>
void init_reference_to_physical_map (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< Node * > &elem_nodes)
virtual void compute_affine_map (const unsigned int dim, const std::vector< Real > &qw, const Elem *elem)
virtual void compute_null_map (const unsigned int dim, const std::vector< Real > &qw)
virtual void compute_map (const unsigned int dim, const std::vector< Real > &qw, const Elem *elem)
virtual void compute_face_map (int dim, const std::vector< Real > &qw, const Elem *side)
void compute_edge_map (int dim, const std::vector< Real > &qw, const Elem *side)
template<unsigned int Dim>
void init_face_shape_functions (const std::vector< Point > &qp, const Elem *side)
template<unsigned int Dim>
void init_edge_shape_functions (const std::vector< Point > &qp, const Elem *edge)
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_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
void print_xyz (std::ostream &os) const
std::vector< std::vector< Real > > & get_psi ()
std::vector< std::vector< Real > > & get_dpsidxi ()
std::vector< std::vector< Real > > & get_dpsideta ()
std::vector< std::vector< Real > > & get_d2psidxi2 ()
std::vector< std::vector< Real > > & get_d2psidxideta ()
std::vector< std::vector< Real > > & get_d2psideta2 ()
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 ()

Static Public Member Functions

static UniquePtr< FEMapbuild (FEType fe_type)

Protected Member Functions

void resize_quadrature_map_vectors (const unsigned int dim, unsigned int n_qp)
Real dxdxi_map (const unsigned int p) const
Real dydxi_map (const unsigned int p) const
Real dzdxi_map (const unsigned int p) const
Real dxdeta_map (const unsigned int p) const
Real dydeta_map (const unsigned int p) const
Real dzdeta_map (const unsigned int p) const
Real dxdzeta_map (const unsigned int p) const
Real dydzeta_map (const unsigned int p) const
Real dzdzeta_map (const unsigned int p) const

Protected Attributes

std::vector< Pointxyz
std::vector< RealGradientdxyzdxi_map
std::vector< RealGradientdxyzdeta_map
std::vector< RealGradientdxyzdzeta_map
std::vector< RealGradientd2xyzdxi2_map
std::vector< RealGradientd2xyzdxideta_map
std::vector< RealGradientd2xyzdeta2_map
std::vector< RealGradientd2xyzdxidzeta_map
std::vector< RealGradientd2xyzdetadzeta_map
std::vector< RealGradientd2xyzdzeta2_map
std::vector< Realdxidx_map
std::vector< Realdxidy_map
std::vector< Realdxidz_map
std::vector< Realdetadx_map
std::vector< Realdetady_map
std::vector< Realdetadz_map
std::vector< Realdzetadx_map
std::vector< Realdzetady_map
std::vector< Realdzetadz_map
std::vector< std::vector< Real > > phi_map
std::vector< std::vector< Real > > dphidxi_map
std::vector< std::vector< Real > > dphideta_map
std::vector< std::vector< Real > > dphidzeta_map
std::vector< std::vector< Real > > d2phidxi2_map
std::vector< std::vector< Real > > d2phidxideta_map
std::vector< std::vector< Real > > d2phidxidzeta_map
std::vector< std::vector< Real > > d2phideta2_map
std::vector< std::vector< Real > > d2phidetadzeta_map
std::vector< std::vector< Real > > d2phidzeta2_map
std::vector< std::vector< Real > > psi_map
std::vector< std::vector< Real > > dpsidxi_map
std::vector< std::vector< Real > > dpsideta_map
std::vector< std::vector< Real > > d2psidxi2_map
std::vector< std::vector< Real > > d2psidxideta_map
std::vector< std::vector< Real > > d2psideta2_map
std::vector< std::vector< Point > > tangents
std::vector< Pointnormals
std::vector< Realcurvatures
std::vector< Realjac
std::vector< RealJxW

Detailed Description

Definition at line 37 of file fe_map.h.


Constructor & Destructor Documentation

Definition at line 38 of file fe_map.C.

Referenced by build().

{}
virtual libMesh::FEMap::~FEMap ( ) [inline, virtual]

Definition at line 42 of file fe_map.h.

{}

Member Function Documentation

UniquePtr< FEMap > libMesh::FEMap::build ( FEType  fe_type) [static]

Definition at line 42 of file fe_map.C.

References libMesh::FEType::family, FEMap(), and libMesh::XYZ.

{
  switch( fe_type.family )
    {
    case XYZ:
      return UniquePtr<FEMap>(new FEXYZMap);

    default:
      return UniquePtr<FEMap>(new FEMap);
    }

  libmesh_error_msg("We'll never get here!");
  return UniquePtr<FEMap>();
}
void libMesh::FEMap::compute_affine_map ( const unsigned int  dim,
const std::vector< Real > &  qw,
const Elem elem 
) [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 690 of file fe_map.C.

References compute_single_point_map(), d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, detadx_map, detady_map, detadz_map, dxidx_map, dxidy_map, dxidz_map, dxyzdeta_map, dxyzdxi_map, dxyzdzeta_map, dzetadx_map, dzetady_map, dzetadz_map, libMesh::Elem::get_node(), jac, JxW, libMesh::libmesh_assert(), libMesh::Elem::n_nodes(), phi_map, resize_quadrature_map_vectors(), libMesh::START_LOG(), and xyz.

Referenced by compute_map().

{
  // Start logging the map computation.
  START_LOG("compute_affine_map()", "FEMap");

  libmesh_assert(elem);

  const unsigned int n_qp = cast_int<unsigned int>(qw.size());

  // Resize the vectors to hold data at the quadrature points
  this->resize_quadrature_map_vectors(dim, n_qp);

  // Determine the nodes contributing to element elem
  std::vector<Node*> elem_nodes(elem->n_nodes(), NULL);
  for (unsigned int i=0; i<elem->n_nodes(); i++)
    elem_nodes[i] = elem->get_node(i);

  // Compute map at quadrature point 0
  this->compute_single_point_map(dim, qw, elem, 0, elem_nodes);

  // Compute xyz at all other quadrature points
  for (unsigned int p=1; p<n_qp; p++)
    {
      xyz[p].zero();
      for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes
        xyz[p].add_scaled        (*elem_nodes[i], phi_map[i][p]    );
    }

  // Copy other map data from quadrature point 0
  for (unsigned int p=1; p<n_qp; p++) // for each extra quadrature point
    {
      dxyzdxi_map[p] = dxyzdxi_map[0];
      dxidx_map[p] = dxidx_map[0];
      dxidy_map[p] = dxidy_map[0];
      dxidz_map[p] = dxidz_map[0];
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      // The map should be affine, so second derivatives are zero
      d2xyzdxi2_map[p] = 0.;
#endif
      if (dim > 1)
        {
          dxyzdeta_map[p] = dxyzdeta_map[0];
          detadx_map[p] = detadx_map[0];
          detady_map[p] = detady_map[0];
          detadz_map[p] = detadz_map[0];
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
          d2xyzdxideta_map[p] = 0.;
          d2xyzdeta2_map[p] = 0.;
#endif
          if (dim > 2)
            {
              dxyzdzeta_map[p] = dxyzdzeta_map[0];
              dzetadx_map[p] = dzetadx_map[0];
              dzetady_map[p] = dzetady_map[0];
              dzetadz_map[p] = dzetadz_map[0];
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
              d2xyzdxidzeta_map[p] = 0.;
              d2xyzdetadzeta_map[p] = 0.;
              d2xyzdzeta2_map[p] = 0.;
#endif
            }
        }
      jac[p] = jac[0];
      JxW[p] = JxW[0] / qw[0] * qw[p];
    }

  STOP_LOG("compute_affine_map()", "FEMap");
}
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 800 of file fe_boundary.C.

References compute_face_map(), curvatures, d2psidxi2_map, d2xyzdeta2_map, d2xyzdxi2_map, d2xyzdxideta_map, dpsidxi_map, dxdxi_map(), dxyzdeta_map, dxyzdxi_map, dydxi_map(), dzdxi_map(), JxW, libMesh::libmesh_assert(), normals, libMesh::Elem::point(), psi_map, libMesh::Real, libMesh::START_LOG(), tangents, and xyz.

{
  libmesh_assert(edge);

  if (dim == 2)
    {
      // A 2D finite element living in either 2D or 3D space.
      // The edges here are the sides of the element, so the
      // (misnamed) compute_face_map function does what we want
      this->compute_face_map(dim, qw, edge);
      return;
    }

  libmesh_assert_equal_to (dim, 3);  // 1D is unnecessary and currently unsupported

  START_LOG("compute_edge_map()", "FEMap");

  // The number of quadrature points.
  const unsigned int n_qp = cast_int<unsigned int>(qw.size());

  // Resize the vectors to hold data at the quadrature points
  this->xyz.resize(n_qp);
  this->dxyzdxi_map.resize(n_qp);
  this->dxyzdeta_map.resize(n_qp);
  this->d2xyzdxi2_map.resize(n_qp);
  this->d2xyzdxideta_map.resize(n_qp);
  this->d2xyzdeta2_map.resize(n_qp);
  this->tangents.resize(n_qp);
  this->normals.resize(n_qp);
  this->curvatures.resize(n_qp);

  this->JxW.resize(n_qp);

  // Clear the entities that will be summed
  for (unsigned int p=0; p<n_qp; p++)
    {
      this->tangents[p].resize(1);
      this->xyz[p].zero();
      this->dxyzdxi_map[p].zero();
      this->dxyzdeta_map[p].zero();
      this->d2xyzdxi2_map[p].zero();
      this->d2xyzdxideta_map[p].zero();
      this->d2xyzdeta2_map[p].zero();
    }

  // compute x, dxdxi at the quadrature points
  for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes
    {
      const Point& edge_point = edge->point(i);

      for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
        {
          this->xyz[p].add_scaled             (edge_point, this->psi_map[i][p]);
          this->dxyzdxi_map[p].add_scaled     (edge_point, this->dpsidxi_map[i][p]);
          this->d2xyzdxi2_map[p].add_scaled   (edge_point, this->d2psidxi2_map[i][p]);
        }
    }

  // Compute the tangents at the quadrature point
  // FIXME: normals (plural!) and curvatures are uncalculated
  for (unsigned int p=0; p<n_qp; p++)
    {
      const Point n  = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]);
      this->tangents[p][0] = this->dxyzdxi_map[p].unit();

      // compute the jacobian at the quadrature points
      const Real the_jac = std::sqrt(this->dxdxi_map(p)*this->dxdxi_map(p) +
                                     this->dydxi_map(p)*this->dydxi_map(p) +
                                     this->dzdxi_map(p)*this->dzdxi_map(p));

      libmesh_assert_greater (the_jac, 0.);

      this->JxW[p] = the_jac*qw[p];
    }

  STOP_LOG("compute_edge_map()", "FEMap");
}
void libMesh::FEMap::compute_face_map ( int  dim,
const std::vector< Real > &  qw,
const Elem side 
) [virtual]

Same as compute_map, but for a side. Useful for boundary integration.

Reimplemented in libMesh::FEXYZMap.

Definition at line 524 of file fe_boundary.C.

References libMesh::TypeVector< T >::cross(), curvatures, d2psideta2_map, d2psidxi2_map, d2psidxideta_map, d2xyzdeta2_map, d2xyzdxi2_map, d2xyzdxideta_map, dpsideta_map, dpsidxi_map, dxdeta_map(), dxdxi_map(), dxyzdeta_map, dxyzdxi_map, dydeta_map(), dydxi_map(), dzdeta_map(), dzdxi_map(), JxW, libMesh::libmesh_assert(), libMesh::Elem::node(), normals, libMesh::Elem::parent(), libMesh::Elem::point(), psi_map, libMesh::Real, libMesh::START_LOG(), tangents, libMesh::TypeVector< T >::unit(), and xyz.

Referenced by compute_edge_map().

{
  libmesh_assert(side);

  START_LOG("compute_face_map()", "FEMap");

  // The number of quadrature points.
  const unsigned int n_qp = cast_int<unsigned int>(qw.size());

  switch (dim)
    {
    case 1:
      {
        // A 1D finite element, currently assumed to be in 1D space
        // This means the boundary is a "0D finite element", a
        // NODEELEM.

        // Resize the vectors to hold data at the quadrature points
        {
          this->xyz.resize(n_qp);
          normals.resize(n_qp);

          this->JxW.resize(n_qp);
        }

        // If we have no quadrature points, there's nothing else to do
        if (!n_qp)
          break;

        // We need to look back at the full edge to figure out the normal
        // vector
        const Elem *elem = side->parent();
        libmesh_assert (elem);
        if (side->node(0) == elem->node(0))
          normals[0] = Point(-1.);
        else
          {
            libmesh_assert_equal_to (side->node(0), elem->node(1));
            normals[0] = Point(1.);
          }

        // Calculate x at the point
        libmesh_assert_equal_to (this->psi_map.size(), 1);
        // In the unlikely event we have multiple quadrature
        // points, they'll be in the same place
        for (unsigned int p=0; p<n_qp; p++)
          {
            this->xyz[p].zero();
            this->xyz[p].add_scaled          (side->point(0), this->psi_map[0][p]);
            normals[p] = normals[0];
            this->JxW[p] = 1.0*qw[p];
          }

        // done computing the map
        break;
      }

    case 2:
      {
        // A 2D finite element living in either 2D or 3D space.
        // This means the boundary is a 1D finite element, i.e.
        // and EDGE2 or EDGE3.
        // Resize the vectors to hold data at the quadrature points
        {
          this->xyz.resize(n_qp);
          this->dxyzdxi_map.resize(n_qp);
          this->d2xyzdxi2_map.resize(n_qp);
          this->tangents.resize(n_qp);
          this->normals.resize(n_qp);
          this->curvatures.resize(n_qp);

          this->JxW.resize(n_qp);
        }

        // Clear the entities that will be summed
        // Compute the tangent & normal at the quadrature point
        for (unsigned int p=0; p<n_qp; p++)
          {
            this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D
            this->xyz[p].zero();
            this->dxyzdxi_map[p].zero();
            this->d2xyzdxi2_map[p].zero();
          }

        // compute x, dxdxi at the quadrature points
        for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes
          {
            const Point& side_point = side->point(i);

            for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
              {
                this->xyz[p].add_scaled          (side_point, this->psi_map[i][p]);
                this->dxyzdxi_map[p].add_scaled  (side_point, this->dpsidxi_map[i][p]);
                this->d2xyzdxi2_map[p].add_scaled(side_point, this->d2psidxi2_map[i][p]);
              }
          }

        // Compute the tangent & normal at the quadrature point
        for (unsigned int p=0; p<n_qp; p++)
          {
            // The first tangent comes from just the edge's Jacobian
            this->tangents[p][0] = this->dxyzdxi_map[p].unit();

#if LIBMESH_DIM == 2
            // For a 2D element living in 2D, the normal is given directly
            // from the entries in the edge Jacobian.
            this->normals[p] = (Point(this->dxyzdxi_map[p](1), -this->dxyzdxi_map[p](0), 0.)).unit();

#elif LIBMESH_DIM == 3
            // For a 2D element living in 3D, there is a second tangent.
            // For the second tangent, we need to refer to the full
            // element's (not just the edge's) Jacobian.
            const Elem *elem = side->parent();
            libmesh_assert(elem);

            // Inverse map xyz[p] to a reference point on the parent...
            Point reference_point = FE<2,LAGRANGE>::inverse_map(elem, this->xyz[p]);

            // Get dxyz/dxi and dxyz/deta from the parent map.
            Point dx_dxi  = FE<2,LAGRANGE>::map_xi (elem, reference_point);
            Point dx_deta = FE<2,LAGRANGE>::map_eta(elem, reference_point);

            // The second tangent vector is formed by crossing these vectors.
            tangents[p][1] = dx_dxi.cross(dx_deta).unit();

            // Finally, the normal in this case is given by crossing these
            // two tangents.
            normals[p] = tangents[p][0].cross(tangents[p][1]).unit();
#endif


            // The curvature is computed via the familiar Frenet formula:
            // curvature = [d^2(x) / d (xi)^2] dot [normal]
            // For a reference, see:
            // F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill, p. 310
            //
            // Note: The sign convention here is different from the
            // 3D case.  Concave-upward curves (smiles) have a positive
            // curvature.  Concave-downward curves (frowns) have a
            // negative curvature.  Be sure to take that into account!
            const Real numerator   = this->d2xyzdxi2_map[p] * this->normals[p];
            const Real denominator = this->dxyzdxi_map[p].size_sq();
            libmesh_assert_not_equal_to (denominator, 0);
            curvatures[p] = numerator / denominator;
          }

        // compute the jacobian at the quadrature points
        for (unsigned int p=0; p<n_qp; p++)
          {
            const Real the_jac = this->dxyzdxi_map[p].size();

            libmesh_assert_greater (the_jac, 0.);

            this->JxW[p] = the_jac*qw[p];
          }

        // done computing the map
        break;
      }



    case 3:
      {
        // A 3D finite element living in 3D space.
        // Resize the vectors to hold data at the quadrature points
        {
          this->xyz.resize(n_qp);
          this->dxyzdxi_map.resize(n_qp);
          this->dxyzdeta_map.resize(n_qp);
          this->d2xyzdxi2_map.resize(n_qp);
          this->d2xyzdxideta_map.resize(n_qp);
          this->d2xyzdeta2_map.resize(n_qp);
          this->tangents.resize(n_qp);
          this->normals.resize(n_qp);
          this->curvatures.resize(n_qp);

          this->JxW.resize(n_qp);
        }

        // Clear the entities that will be summed
        for (unsigned int p=0; p<n_qp; p++)
          {
            this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D
            this->xyz[p].zero();
            this->dxyzdxi_map[p].zero();
            this->dxyzdeta_map[p].zero();
            this->d2xyzdxi2_map[p].zero();
            this->d2xyzdxideta_map[p].zero();
            this->d2xyzdeta2_map[p].zero();
          }

        // compute x, dxdxi at the quadrature points
        for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes
          {
            const Point& side_point = side->point(i);

            for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
              {
                this->xyz[p].add_scaled         (side_point, this->psi_map[i][p]);
                this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]);
                this->dxyzdeta_map[p].add_scaled(side_point, this->dpsideta_map[i][p]);
                this->d2xyzdxi2_map[p].add_scaled   (side_point, this->d2psidxi2_map[i][p]);
                this->d2xyzdxideta_map[p].add_scaled(side_point, this->d2psidxideta_map[i][p]);
                this->d2xyzdeta2_map[p].add_scaled  (side_point, this->d2psideta2_map[i][p]);
              }
          }

        // Compute the tangents, normal, and curvature at the quadrature point
        for (unsigned int p=0; p<n_qp; p++)
          {
            const Point n  = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]);
            this->normals[p]     = n.unit();
            this->tangents[p][0] = this->dxyzdxi_map[p].unit();
            this->tangents[p][1] = n.cross(this->dxyzdxi_map[p]).unit();

            // Compute curvature using the typical nomenclature
            // of the first and second fundamental forms.
            // For reference, see:
            // 1) http://mathworld.wolfram.com/MeanCurvature.html
            //    (note -- they are using inward normal)
            // 2) F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill
            const Real L  = -this->d2xyzdxi2_map[p]    * this->normals[p];
            const Real M  = -this->d2xyzdxideta_map[p] * this->normals[p];
            const Real N  = -this->d2xyzdeta2_map[p]   * this->normals[p];
            const Real E  =  this->dxyzdxi_map[p].size_sq();
            const Real F  =  this->dxyzdxi_map[p]      * this->dxyzdeta_map[p];
            const Real G  =  this->dxyzdeta_map[p].size_sq();

            const Real numerator   = E*N -2.*F*M + G*L;
            const Real denominator = E*G - F*F;
            libmesh_assert_not_equal_to (denominator, 0.);
            curvatures[p] = 0.5*numerator/denominator;
          }

        // compute the jacobian at the quadrature points, see
        // http://sp81.msi.umn.edu:999/fluent/fidap/help/theory/thtoc.htm
        for (unsigned int p=0; p<n_qp; p++)
          {
            const Real g11 = (dxdxi_map(p)*dxdxi_map(p) +
                              dydxi_map(p)*dydxi_map(p) +
                              dzdxi_map(p)*dzdxi_map(p));

            const Real g12 = (dxdxi_map(p)*dxdeta_map(p) +
                              dydxi_map(p)*dydeta_map(p) +
                              dzdxi_map(p)*dzdeta_map(p));

            const Real g21 = g12;

            const Real g22 = (dxdeta_map(p)*dxdeta_map(p) +
                              dydeta_map(p)*dydeta_map(p) +
                              dzdeta_map(p)*dzdeta_map(p));


            const Real the_jac = std::sqrt(g11*g22 - g12*g21);

            libmesh_assert_greater (the_jac, 0.);

            this->JxW[p] = the_jac*qw[p];
          }

        // done computing the map
        break;
      }


    default:
      libmesh_error_msg("Invalid dimension dim = " << dim);
    }
  STOP_LOG("compute_face_map()", "FEMap");
}
void libMesh::FEMap::compute_map ( const unsigned int  dim,
const std::vector< Real > &  qw,
const Elem elem 
) [virtual]

Compute the jacobian and some other additional data fields. Takes the integration weights as input, along with a pointer to the element.

Definition at line 818 of file fe_map.C.

References compute_affine_map(), compute_null_map(), compute_single_point_map(), libMesh::MeshTools::Subdivision::find_one_ring(), libMesh::Elem::get_node(), libMesh::Elem::has_affine_map(), libMesh::libmesh_assert(), libMesh::Elem::n_nodes(), resize_quadrature_map_vectors(), libMesh::START_LOG(), libMesh::TRI3SUBDIVISION, and libMesh::Elem::type().

{
  if (!elem)
    {
      compute_null_map(dim, qw);
      return;
    }

  if (elem->has_affine_map())
    {
      compute_affine_map(dim, qw, elem);
      return;
    }

  // Start logging the map computation.
  START_LOG("compute_map()", "FEMap");

  libmesh_assert(elem);

  const unsigned int n_qp = cast_int<unsigned int>(qw.size());

  // Resize the vectors to hold data at the quadrature points
  this->resize_quadrature_map_vectors(dim, n_qp);

  // Determine the nodes contributing to element elem
  std::vector<Node*> elem_nodes;
  if (elem->type() == TRI3SUBDIVISION)
    {
      // Subdivision surface FE require the 1-ring around elem
      libmesh_assert_equal_to (dim, 2);
      const Tri3Subdivision* sd_elem = static_cast<const Tri3Subdivision*>(elem);
      MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
    }
  else
    {
      // All other FE use only the nodes of elem itself
      elem_nodes.resize(elem->n_nodes(), NULL);
      for (unsigned int i=0; i<elem->n_nodes(); i++)
        elem_nodes[i] = elem->get_node(i);
    }

  // Compute map at all quadrature points
  for (unsigned int p=0; p!=n_qp; p++)
    this->compute_single_point_map(dim, qw, elem, p, elem_nodes);

  // Stop logging the map computation.
  STOP_LOG("compute_map()", "FEMap");
}
void libMesh::FEMap::compute_null_map ( const unsigned int  dim,
const std::vector< Real > &  qw 
) [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 763 of file fe_map.C.

References d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, detadx_map, detady_map, detadz_map, dxidx_map, dxidy_map, dxidz_map, dxyzdeta_map, dxyzdxi_map, dxyzdzeta_map, dzetadx_map, dzetady_map, dzetadz_map, jac, JxW, resize_quadrature_map_vectors(), libMesh::START_LOG(), and xyz.

Referenced by compute_map().

{
  // Start logging the map computation.
  START_LOG("compute_null_map()", "FEMap");

  const unsigned int n_qp = cast_int<unsigned int>(qw.size());

  // Resize the vectors to hold data at the quadrature points
  this->resize_quadrature_map_vectors(dim, n_qp);

  // Compute "fake" xyz
  for (unsigned int p=1; p<n_qp; p++)
    {
      xyz[p].zero();

      dxyzdxi_map[p] = 0;
      dxidx_map[p] = 0;
      dxidy_map[p] = 0;
      dxidz_map[p] = 0;
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      d2xyzdxi2_map[p] = 0;
#endif
      if (dim > 1)
        {
          dxyzdeta_map[p] = 0;
          detadx_map[p] = 0;
          detady_map[p] = 0;
          detadz_map[p] = 0;
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
          d2xyzdxideta_map[p] = 0.;
          d2xyzdeta2_map[p] = 0.;
#endif
          if (dim > 2)
            {
              dxyzdzeta_map[p] = 0;
              dzetadx_map[p] = 0;
              dzetady_map[p] = 0;
              dzetadz_map[p] = 0;
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
              d2xyzdxidzeta_map[p] = 0;
              d2xyzdetadzeta_map[p] = 0;
              d2xyzdzeta2_map[p] = 0;
#endif
            }
        }
      jac[p] = 1;
      JxW[p] = qw[p];
    }

  STOP_LOG("compute_null_map()", "FEMap");
}
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< Node * > &  elem_nodes 
)

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.

Definition at line 312 of file fe_map.C.

References d2phideta2_map, d2phidetadzeta_map, d2phidxi2_map, d2phidxideta_map, d2phidxidzeta_map, d2phidzeta2_map, d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, detadx_map, detady_map, detadz_map, 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::DofObject::id(), jac, JxW, libMesh::libmesh_assert(), phi_map, libMesh::Real, and xyz.

Referenced by compute_affine_map(), and compute_map().

{
  libmesh_assert(elem);
  libmesh_assert_equal_to(phi_map.size(), elem_nodes.size());

  switch (dim)
    {
      //--------------------------------------------------------------------
      // 0D
    case 0:
      {
        libmesh_assert(elem_nodes[0]);
        xyz[p] = *elem_nodes[0];
        jac[p] = 1.0;
        JxW[p] = qw[p];
        break;
      }

      //--------------------------------------------------------------------
      // 1D
    case 1:
      {
        // Clear the entities that will be summed
        xyz[p].zero();
        dxyzdxi_map[p].zero();
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
        d2xyzdxi2_map[p].zero();
#endif

        // compute x, dx, d2x at the quadrature point
        for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes
          {
            // Reference to the point, helps eliminate
            // exessive temporaries in the inner loop
            libmesh_assert(elem_nodes[i]);
            const Point& elem_point = *elem_nodes[i];

            xyz[p].add_scaled          (elem_point, phi_map[i][p]    );
            dxyzdxi_map[p].add_scaled  (elem_point, dphidxi_map[i][p]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
            d2xyzdxi2_map[p].add_scaled(elem_point, d2phidxi2_map[i][p]);
#endif
          }

        // Compute the jacobian
        //
        // 1D elements can live in 2D or 3D space.
        // The transformation matrix from local->global
        // coordinates is
        //
        // T = | dx/dxi |
        //     | dy/dxi |
        //     | dz/dxi |
        //
        // The generalized determinant of T (from the
        // so-called "normal" eqns.) is
        // jac = "det(T)" = sqrt(det(T'T))
        //
        // where T'= transpose of T, so
        //
        // jac = sqrt( (dx/dxi)^2 + (dy/dxi)^2 + (dz/dxi)^2 )
        jac[p] = dxyzdxi_map[p].size();

        if (jac[p] <= 0.)
          libmesh_error_msg("ERROR: negative Jacobian: " << jac[p] << " in element " << elem->id());

        // The inverse Jacobian entries also come from the
        // generalized inverse of T (see also the 2D element
        // living in 3D code).
        const Real jacm2 = 1./jac[p]/jac[p];
        dxidx_map[p] = jacm2*dxdxi_map(p);
#if LIBMESH_DIM > 1
        dxidy_map[p] = jacm2*dydxi_map(p);
#endif
#if LIBMESH_DIM > 2
        dxidz_map[p] = jacm2*dzdxi_map(p);
#endif

        JxW[p] = jac[p]*qw[p];

        // done computing the map
        break;
      }


      //--------------------------------------------------------------------
      // 2D
    case 2:
      {
        //------------------------------------------------------------------
        // Compute the (x,y) values at the quadrature points,
        // the Jacobian at the quadrature points

        xyz[p].zero();

        dxyzdxi_map[p].zero();
        dxyzdeta_map[p].zero();
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
        d2xyzdxi2_map[p].zero();
        d2xyzdxideta_map[p].zero();
        d2xyzdeta2_map[p].zero();
#endif


        // compute (x,y) at the quadrature points, derivatives once
        for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes
          {
            // Reference to the point, helps eliminate
            // exessive temporaries in the inner loop
            libmesh_assert(elem_nodes[i]);
            const Point& elem_point = *elem_nodes[i];

            xyz[p].add_scaled          (elem_point, phi_map[i][p]     );

            dxyzdxi_map[p].add_scaled      (elem_point, dphidxi_map[i][p] );
            dxyzdeta_map[p].add_scaled     (elem_point, dphideta_map[i][p]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
            d2xyzdxi2_map[p].add_scaled    (elem_point, d2phidxi2_map[i][p]);
            d2xyzdxideta_map[p].add_scaled (elem_point, d2phidxideta_map[i][p]);
            d2xyzdeta2_map[p].add_scaled   (elem_point, d2phideta2_map[i][p]);
#endif
          }

        // compute the jacobian once
        const Real dx_dxi = dxdxi_map(p),
          dx_deta = dxdeta_map(p),
          dy_dxi = dydxi_map(p),
          dy_deta = dydeta_map(p);

#if LIBMESH_DIM == 2
        // Compute the Jacobian.  This assumes the 2D face
        // lives in 2D space
        //
        // Symbolically, the matrix determinant is
        //
        //         | dx/dxi  dx/deta |
        // jac =   | dy/dxi  dy/deta |
        //
        // jac = dx/dxi*dy/deta - dx/deta*dy/dxi
        jac[p] = (dx_dxi*dy_deta - dx_deta*dy_dxi);

        if (jac[p] <= 0.)
          libmesh_error_msg("ERROR: negative Jacobian: " << jac[p] << " in element " << elem->id());

        JxW[p] = jac[p]*qw[p];

        // Compute the shape function derivatives wrt x,y at the
        // quadrature points
        const Real inv_jac = 1./jac[p];

        dxidx_map[p]  =  dy_deta*inv_jac; //dxi/dx  =  (1/J)*dy/deta
        dxidy_map[p]  = -dx_deta*inv_jac; //dxi/dy  = -(1/J)*dx/deta
        detadx_map[p] = -dy_dxi* inv_jac; //deta/dx = -(1/J)*dy/dxi
        detady_map[p] =  dx_dxi* inv_jac; //deta/dy =  (1/J)*dx/dxi

        dxidz_map[p] = detadz_map[p] = 0.;
#else

        const Real dz_dxi = dzdxi_map(p),
          dz_deta = dzdeta_map(p);

        // Compute the Jacobian.  This assumes a 2D face in
        // 3D space.
        //
        // The transformation matrix T from local to global
        // coordinates is
        //
        //         | dx/dxi  dx/deta |
        //     T = | dy/dxi  dy/deta |
        //         | dz/dxi  dz/deta |
        // note det(T' T) = det(T')det(T) = det(T)det(T)
        // so det(T) = std::sqrt(det(T' T))
        //
        //----------------------------------------------
        // Notes:
        //
        //       dX = R dXi -> R'dX = R'R dXi
        // (R^-1)dX =   dXi    [(R'R)^-1 R']dX = dXi
        //
        // so R^-1 = (R'R)^-1 R'
        //
        // and R^-1 R = (R'R)^-1 R'R = I.
        //
        const Real g11 = (dx_dxi*dx_dxi +
                          dy_dxi*dy_dxi +
                          dz_dxi*dz_dxi);

        const Real g12 = (dx_dxi*dx_deta +
                          dy_dxi*dy_deta +
                          dz_dxi*dz_deta);

        const Real g21 = g12;

        const Real g22 = (dx_deta*dx_deta +
                          dy_deta*dy_deta +
                          dz_deta*dz_deta);

        const Real det = (g11*g22 - g12*g21);

        if (det <= 0.)
          libmesh_error_msg("ERROR: negative Jacobian! " << " in element " << elem->id());

        const Real inv_det = 1./det;
        jac[p] = std::sqrt(det);

        JxW[p] = jac[p]*qw[p];

        const Real g11inv =  g22*inv_det;
        const Real g12inv = -g12*inv_det;
        const Real g21inv = -g21*inv_det;
        const Real g22inv =  g11*inv_det;

        dxidx_map[p]  = g11inv*dx_dxi + g12inv*dx_deta;
        dxidy_map[p]  = g11inv*dy_dxi + g12inv*dy_deta;
        dxidz_map[p]  = g11inv*dz_dxi + g12inv*dz_deta;

        detadx_map[p] = g21inv*dx_dxi + g22inv*dx_deta;
        detady_map[p] = g21inv*dy_dxi + g22inv*dy_deta;
        detadz_map[p] = g21inv*dz_dxi + g22inv*dz_deta;

#endif
        // done computing the map
        break;
      }



      //--------------------------------------------------------------------
      // 3D
    case 3:
      {
        //------------------------------------------------------------------
        // Compute the (x,y,z) values at the quadrature points,
        // the Jacobian at the quadrature point

        // Clear the entities that will be summed
        xyz[p].zero           ();
        dxyzdxi_map[p].zero   ();
        dxyzdeta_map[p].zero  ();
        dxyzdzeta_map[p].zero ();
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
        d2xyzdxi2_map[p].zero();
        d2xyzdxideta_map[p].zero();
        d2xyzdxidzeta_map[p].zero();
        d2xyzdeta2_map[p].zero();
        d2xyzdetadzeta_map[p].zero();
        d2xyzdzeta2_map[p].zero();
#endif


        // compute (x,y,z) at the quadrature points,
        // dxdxi,   dydxi,   dzdxi,
        // dxdeta,  dydeta,  dzdeta,
        // dxdzeta, dydzeta, dzdzeta  all once
        for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes
          {
            // Reference to the point, helps eliminate
            // exessive temporaries in the inner loop
            libmesh_assert(elem_nodes[i]);
            const Point& elem_point = *elem_nodes[i];

            xyz[p].add_scaled           (elem_point, phi_map[i][p]      );
            dxyzdxi_map[p].add_scaled   (elem_point, dphidxi_map[i][p]  );
            dxyzdeta_map[p].add_scaled  (elem_point, dphideta_map[i][p] );
            dxyzdzeta_map[p].add_scaled (elem_point, dphidzeta_map[i][p]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
            d2xyzdxi2_map[p].add_scaled      (elem_point,
                                              d2phidxi2_map[i][p]);
            d2xyzdxideta_map[p].add_scaled   (elem_point,
                                              d2phidxideta_map[i][p]);
            d2xyzdxidzeta_map[p].add_scaled  (elem_point,
                                              d2phidxidzeta_map[i][p]);
            d2xyzdeta2_map[p].add_scaled     (elem_point,
                                              d2phideta2_map[i][p]);
            d2xyzdetadzeta_map[p].add_scaled (elem_point,
                                              d2phidetadzeta_map[i][p]);
            d2xyzdzeta2_map[p].add_scaled    (elem_point,
                                              d2phidzeta2_map[i][p]);
#endif
          }

        // compute the jacobian
        const Real
          dx_dxi   = dxdxi_map(p),   dy_dxi   = dydxi_map(p),   dz_dxi   = dzdxi_map(p),
          dx_deta  = dxdeta_map(p),  dy_deta  = dydeta_map(p),  dz_deta  = dzdeta_map(p),
          dx_dzeta = dxdzeta_map(p), dy_dzeta = dydzeta_map(p), dz_dzeta = dzdzeta_map(p);

        // Symbolically, the matrix determinant is
        //
        //         | dx/dxi   dy/dxi   dz/dxi   |
        // jac =   | dx/deta  dy/deta  dz/deta  |
        //         | dx/dzeta dy/dzeta dz/dzeta |
        //
        // jac = dx/dxi*(dy/deta*dz/dzeta - dz/deta*dy/dzeta) +
        //       dy/dxi*(dz/deta*dx/dzeta - dx/deta*dz/dzeta) +
        //       dz/dxi*(dx/deta*dy/dzeta - dy/deta*dx/dzeta)

        jac[p] = (dx_dxi*(dy_deta*dz_dzeta - dz_deta*dy_dzeta)  +
                  dy_dxi*(dz_deta*dx_dzeta - dx_deta*dz_dzeta)  +
                  dz_dxi*(dx_deta*dy_dzeta - dy_deta*dx_dzeta));

        if (jac[p] <= 0.)
          libmesh_error_msg("ERROR: negative Jacobian: " << jac[p] << " in element " << elem->id());

        JxW[p] = jac[p]*qw[p];

        // Compute the shape function derivatives wrt x,y at the
        // quadrature points
        const Real inv_jac  = 1./jac[p];

        dxidx_map[p]   = (dy_deta*dz_dzeta - dz_deta*dy_dzeta)*inv_jac;
        dxidy_map[p]   = (dz_deta*dx_dzeta - dx_deta*dz_dzeta)*inv_jac;
        dxidz_map[p]   = (dx_deta*dy_dzeta - dy_deta*dx_dzeta)*inv_jac;

        detadx_map[p]  = (dz_dxi*dy_dzeta  - dy_dxi*dz_dzeta )*inv_jac;
        detady_map[p]  = (dx_dxi*dz_dzeta  - dz_dxi*dx_dzeta )*inv_jac;
        detadz_map[p]  = (dy_dxi*dx_dzeta  - dx_dxi*dy_dzeta )*inv_jac;

        dzetadx_map[p] = (dy_dxi*dz_deta   - dz_dxi*dy_deta  )*inv_jac;
        dzetady_map[p] = (dz_dxi*dx_deta   - dx_dxi*dz_deta  )*inv_jac;
        dzetadz_map[p] = (dx_dxi*dy_deta   - dy_dxi*dx_deta  )*inv_jac;

        // done computing the map
        break;
      }

    default:
      libmesh_error_msg("Invalid dim = " << dim);
    }
}
Real libMesh::FEMap::dxdeta_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the x value of the pth entry of the dxzydeta_map.

Definition at line 472 of file fe_map.h.

References dxyzdeta_map.

Referenced by libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and compute_single_point_map().

{ return dxyzdeta_map[p](0); }
Real libMesh::FEMap::dxdxi_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the x value of the pth entry of the dxzydxi_map.

Definition at line 451 of file fe_map.h.

References dxyzdxi_map.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and compute_single_point_map().

{ return dxyzdxi_map[p](0); }
Real libMesh::FEMap::dxdzeta_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the x value of the pth entry of the dxzydzeta_map.

Definition at line 493 of file fe_map.h.

References dxyzdzeta_map.

Referenced by compute_single_point_map().

{ return dxyzdzeta_map[p](0); }
Real libMesh::FEMap::dydeta_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the y value of the pth entry of the dxzydeta_map.

Definition at line 479 of file fe_map.h.

References dxyzdeta_map.

Referenced by libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and compute_single_point_map().

{ return dxyzdeta_map[p](1); }
Real libMesh::FEMap::dydxi_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the y value of the pth entry of the dxzydxi_map.

Definition at line 458 of file fe_map.h.

References dxyzdxi_map.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and compute_single_point_map().

{ return dxyzdxi_map[p](1); }
Real libMesh::FEMap::dydzeta_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the y value of the pth entry of the dxzydzeta_map.

Definition at line 500 of file fe_map.h.

References dxyzdzeta_map.

Referenced by compute_single_point_map().

{ return dxyzdzeta_map[p](1); }
Real libMesh::FEMap::dzdeta_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the z value of the pth entry of the dxzydeta_map.

Definition at line 486 of file fe_map.h.

References dxyzdeta_map.

Referenced by libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and compute_single_point_map().

{ return dxyzdeta_map[p](2); }
Real libMesh::FEMap::dzdxi_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the z value of the pth entry of the dxzydxi_map.

Definition at line 465 of file fe_map.h.

References dxyzdxi_map.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and compute_single_point_map().

{ return dxyzdxi_map[p](2); }
Real libMesh::FEMap::dzdzeta_map ( const unsigned int  p) const [inline, protected]

Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the z value of the pth entry of the dxzydzeta_map.

Definition at line 507 of file fe_map.h.

References dxyzdzeta_map.

Referenced by compute_single_point_map().

{ return dxyzdzeta_map[p](2); }
const std::vector<Real>& libMesh::FEMap::get_curvatures ( ) const [inline]
Returns:
the curvatures for use in face integration.

Definition at line 314 of file fe_map.h.

References curvatures.

  { return curvatures;}
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phideta2_map ( ) [inline]
Returns:
the reference to physical map 2nd derivative

Definition at line 414 of file fe_map.h.

References d2phideta2_map.

  { return d2phideta2_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidetadzeta_map ( ) [inline]
Returns:
the reference to physical map 2nd derivative

Definition at line 420 of file fe_map.h.

References d2phidetadzeta_map.

  { return d2phidetadzeta_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidxi2_map ( ) [inline]
Returns:
the reference to physical map 2nd derivative

Definition at line 396 of file fe_map.h.

References d2phidxi2_map.

  { return d2phidxi2_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidxideta_map ( ) [inline]
Returns:
the reference to physical map 2nd derivative

Definition at line 402 of file fe_map.h.

References d2phidxideta_map.

  { return d2phidxideta_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidxidzeta_map ( ) [inline]
Returns:
the reference to physical map 2nd derivative

Definition at line 408 of file fe_map.h.

References d2phidxidzeta_map.

  { return d2phidxidzeta_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidzeta2_map ( ) [inline]
Returns:
the reference to physical map 2nd derivative

Definition at line 426 of file fe_map.h.

References d2phidzeta2_map.

  { return d2phidzeta2_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2psideta2 ( ) [inline]
Returns:
the reference to physical map 2nd derivative for the side/edge

Definition at line 365 of file fe_map.h.

References d2psideta2_map.

  { return d2psideta2_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2psidxi2 ( ) [inline]
Returns:
the reference to physical map 2nd derivative for the side/edge

Definition at line 353 of file fe_map.h.

References d2psidxi2_map.

  { return d2psidxi2_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_d2psidxideta ( ) [inline]
Returns:
the reference to physical map 2nd derivative for the side/edge

Definition at line 359 of file fe_map.h.

References d2psidxideta_map.

  { return d2psidxideta_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdeta2 ( ) const [inline]
Returns:
the second partial derivatives in eta.

Definition at line 171 of file fe_map.h.

References d2xyzdeta2_map.

  { return d2xyzdeta2_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdetadzeta ( ) const [inline]
Returns:
the second partial derivatives in eta-zeta.

Definition at line 201 of file fe_map.h.

References d2xyzdetadzeta_map.

  { return d2xyzdetadzeta_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdxi2 ( ) const [inline]
Returns:
the second partial derivatives in xi.

Definition at line 165 of file fe_map.h.

References d2xyzdxi2_map.

  { return d2xyzdxi2_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdxideta ( ) const [inline]
Returns:
the second partial derivatives in xi-eta.

Definition at line 187 of file fe_map.h.

References d2xyzdxideta_map.

  { return d2xyzdxideta_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdxidzeta ( ) const [inline]
Returns:
the second partial derivatives in xi-zeta.

Definition at line 195 of file fe_map.h.

References d2xyzdxidzeta_map.

  { return d2xyzdxidzeta_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdzeta2 ( ) const [inline]
Returns:
the second partial derivatives in zeta.

Definition at line 179 of file fe_map.h.

References d2xyzdzeta2_map.

  { return d2xyzdzeta2_map; }
const std::vector<std::vector<Real> >& libMesh::FEMap::get_dphideta_map ( ) const [inline]
Returns:
the reference to physical map derivative

Definition at line 290 of file fe_map.h.

References dphideta_map.

  { return dphideta_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_dphideta_map ( ) [inline]
Returns:
the reference to physical map derivative

Definition at line 383 of file fe_map.h.

References dphideta_map.

  { return dphideta_map; }
const std::vector<std::vector<Real> >& libMesh::FEMap::get_dphidxi_map ( ) const [inline]
Returns:
the reference to physical map derivative

Definition at line 284 of file fe_map.h.

References dphidxi_map.

  { return dphidxi_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_dphidxi_map ( ) [inline]
Returns:
the reference to physical map derivative

Definition at line 377 of file fe_map.h.

References dphidxi_map.

  { return dphidxi_map; }
const std::vector<std::vector<Real> >& libMesh::FEMap::get_dphidzeta_map ( ) const [inline]
Returns:
the reference to physical map derivative

Definition at line 296 of file fe_map.h.

References dphidzeta_map.

  { return dphidzeta_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_dphidzeta_map ( ) [inline]
Returns:
the reference to physical map derivative

Definition at line 389 of file fe_map.h.

References dphidzeta_map.

  { return dphidzeta_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_dpsideta ( ) [inline]
Returns:
the reference to physical map derivative for the side/edge

Definition at line 347 of file fe_map.h.

References dpsideta_map.

  { return dpsideta_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_dpsidxi ( ) [inline]
Returns:
the reference to physical map derivative for the side/edge

Definition at line 341 of file fe_map.h.

References dpsidxi_map.

  { return dpsidxi_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_dxyzdeta ( ) const [inline]
Returns:
the element tangents in eta-direction at the quadrature points.

Definition at line 152 of file fe_map.h.

References dxyzdeta_map.

Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl().

  { return dxyzdeta_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_dxyzdxi ( ) const [inline]
Returns:
the element tangents in xi-direction at the quadrature points.

Definition at line 145 of file fe_map.h.

References dxyzdxi_map.

Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl().

  { return dxyzdxi_map; }
const std::vector<RealGradient>& libMesh::FEMap::get_dxyzdzeta ( ) const [inline]
Returns:
the element tangents in zeta-direction at the quadrature points.

Definition at line 159 of file fe_map.h.

References dxyzdzeta_map.

Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl().

  { return dxyzdzeta_map; }
const std::vector<Real>& libMesh::FEMap::get_jacobian ( ) const [inline]
Returns:
the element Jacobian for each quadrature point.

Definition at line 131 of file fe_map.h.

References jac.

Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl().

  { return jac; }
const std::vector<Real>& libMesh::FEMap::get_JxW ( ) const [inline]
Returns:
the element Jacobian times the quadrature weight for each quadrature point.

Definition at line 138 of file fe_map.h.

References JxW.

  { return JxW; }
std::vector<Real>& libMesh::FEMap::get_JxW ( ) [inline]
Returns:
writable reference to the element Jacobian times the quadrature weight for each quadrature point.

Definition at line 436 of file fe_map.h.

References JxW.

  { return JxW; }
const std::vector<Point>& libMesh::FEMap::get_normals ( ) const [inline]
Returns:
the normal vectors for face integration.

Definition at line 308 of file fe_map.h.

References normals.

  { return normals; }
const std::vector<std::vector<Real> >& libMesh::FEMap::get_phi_map ( ) const [inline]
Returns:
the reference to physical map for the element

Definition at line 278 of file fe_map.h.

References phi_map.

  { return phi_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_phi_map ( ) [inline]
Returns:
the reference to physical map for the element

Definition at line 371 of file fe_map.h.

References phi_map.

  { return phi_map; }
const std::vector<std::vector<Real> >& libMesh::FEMap::get_psi ( ) const [inline]
Returns:
the reference to physical map for the side/edge

Definition at line 272 of file fe_map.h.

References psi_map.

  { return psi_map; }
std::vector<std::vector<Real> >& libMesh::FEMap::get_psi ( ) [inline]
Returns:
the reference to physical map for the side/edge

Definition at line 335 of file fe_map.h.

References psi_map.

  { return psi_map; }
const std::vector<std::vector<Point> >& libMesh::FEMap::get_tangents ( ) const [inline]
Returns:
the tangent vectors for face integration.

Definition at line 302 of file fe_map.h.

References tangents.

  { return tangents; }
const std::vector<Point>& libMesh::FEMap::get_xyz ( ) const [inline]
Returns:
the xyz spatial locations of the quadrature points on the element.

Definition at line 125 of file fe_map.h.

References xyz.

  { return xyz; }
template<unsigned int Dim>
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.

Start logging the shape function initialization

Stop logging the shape function initialization

Definition at line 469 of file fe_boundary.C.

References d2psidxi2_map, libMesh::Elem::default_order(), dpsidxi_map, libMesh::libmesh_assert(), psi_map, libMesh::START_LOG(), and libMesh::Elem::type().

{
  libmesh_assert(edge);

  START_LOG("init_edge_shape_functions()", "FEMap");

  // The element type and order to use in
  // the map
  const Order    mapping_order     (edge->default_order());
  const ElemType mapping_elem_type (edge->type());

  // The number of quadrature points.
  const unsigned int n_qp = cast_int<unsigned int>(qp.size());

  const unsigned int n_mapping_shape_functions =
    FE<Dim,LAGRANGE>::n_shape_functions (mapping_elem_type,
                                         mapping_order);

  // resize the vectors to hold current data
  // Psi are the shape functions used for the FE mapping
  this->psi_map.resize        (n_mapping_shape_functions);
  this->dpsidxi_map.resize    (n_mapping_shape_functions);
  this->d2psidxi2_map.resize  (n_mapping_shape_functions);

  for (unsigned int i=0; i<n_mapping_shape_functions; i++)
    {
      // Allocate space to store the values of the shape functions
      // and their first and second derivatives at the quadrature points.
      this->psi_map[i].resize        (n_qp);
      this->dpsidxi_map[i].resize    (n_qp);
      this->d2psidxi2_map[i].resize  (n_qp);

      // Compute the value of shape function i, and its first and
      // second derivatives at quadrature point p
      // (Lagrange shape functions are used for the mapping)
      for (unsigned int p=0; p<n_qp; p++)
        {
          this->psi_map[i][p]        = FE<1,LAGRANGE>::shape             (mapping_elem_type, mapping_order, i,    qp[p]);
          this->dpsidxi_map[i][p]    = FE<1,LAGRANGE>::shape_deriv       (mapping_elem_type, mapping_order, i, 0, qp[p]);
          this->d2psidxi2_map[i][p]  = FE<1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 0, qp[p]);
        }
    }

  STOP_LOG("init_edge_shape_functions()", "FEMap");
}
template<unsigned int Dim>
void libMesh::FEMap::init_face_shape_functions ( const std::vector< Point > &  qp,
const Elem side 
)

Initalizes the reference to physical element map for a side. This is used for boundary integration.

Start logging the shape function initialization

Stop logging the shape function initialization

Definition at line 380 of file fe_boundary.C.

References d2psideta2_map, d2psidxi2_map, d2psidxideta_map, libMesh::Elem::default_order(), dpsideta_map, dpsidxi_map, libMesh::libmesh_assert(), psi_map, libMesh::START_LOG(), and libMesh::Elem::type().

{
  libmesh_assert(side);

  START_LOG("init_face_shape_functions()", "FEMap");

  // The element type and order to use in
  // the map
  const Order    mapping_order     (side->default_order());
  const ElemType mapping_elem_type (side->type());

  // The number of quadrature points.
  const unsigned int n_qp = cast_int<unsigned int>(qp.size());

  const unsigned int n_mapping_shape_functions =
    FE<Dim,LAGRANGE>::n_shape_functions (mapping_elem_type,
                                         mapping_order);

  // resize the vectors to hold current data
  // Psi are the shape functions used for the FE mapping
  this->psi_map.resize        (n_mapping_shape_functions);

  if (Dim > 1)
    {
      this->dpsidxi_map.resize    (n_mapping_shape_functions);
      this->d2psidxi2_map.resize  (n_mapping_shape_functions);
    }

  if (Dim == 3)
    {
      this->dpsideta_map.resize     (n_mapping_shape_functions);
      this->d2psidxideta_map.resize (n_mapping_shape_functions);
      this->d2psideta2_map.resize   (n_mapping_shape_functions);
    }

  for (unsigned int i=0; i<n_mapping_shape_functions; i++)
    {
      // Allocate space to store the values of the shape functions
      // and their first and second derivatives at the quadrature points.
      this->psi_map[i].resize        (n_qp);
      if (Dim > 1)
        {
          this->dpsidxi_map[i].resize    (n_qp);
          this->d2psidxi2_map[i].resize  (n_qp);
        }
      if (Dim == 3)
        {
          this->dpsideta_map[i].resize     (n_qp);
          this->d2psidxideta_map[i].resize (n_qp);
          this->d2psideta2_map[i].resize   (n_qp);
        }

      // Compute the value of shape function i, and its first and
      // second derivatives at quadrature point p
      // (Lagrange shape functions are used for the mapping)
      for (unsigned int p=0; p<n_qp; p++)
        {
          this->psi_map[i][p]        = FE<Dim-1,LAGRANGE>::shape             (mapping_elem_type, mapping_order, i,    qp[p]);
          if (Dim > 1)
            {
              this->dpsidxi_map[i][p]    = FE<Dim-1,LAGRANGE>::shape_deriv       (mapping_elem_type, mapping_order, i, 0, qp[p]);
              this->d2psidxi2_map[i][p]  = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 0, qp[p]);
            }
          // libMesh::out << "this->d2psidxi2_map["<<i<<"][p]=" << d2psidxi2_map[i][p] << std::endl;

          // If we are in 3D, then our sides are 2D faces.
          // For the second derivatives, we must also compute the cross
          // derivative d^2() / dxi deta
          if (Dim == 3)
            {
              this->dpsideta_map[i][p]     = FE<Dim-1,LAGRANGE>::shape_deriv       (mapping_elem_type, mapping_order, i, 1, qp[p]);
              this->d2psidxideta_map[i][p] = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 1, qp[p]);
              this->d2psideta2_map[i][p]   = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 2, qp[p]);
            }
        }
    }


  STOP_LOG("init_face_shape_functions()", "FEMap");
}
template<unsigned int Dim>
template void libMesh::FEMap::init_reference_to_physical_map< 3 > ( const std::vector< Point > &  qp,
const Elem elem 
)

Definition at line 60 of file fe_map.C.

References d2phideta2_map, d2phidetadzeta_map, d2phidxi2_map, d2phidxideta_map, d2phidxidzeta_map, d2phidzeta2_map, libMesh::Elem::default_order(), dphideta_map, dphidxi_map, dphidzeta_map, libMesh::Elem::is_linear(), phi_map, libMesh::START_LOG(), and libMesh::Elem::type().

{
  // Start logging the reference->physical map initialization
  START_LOG("init_reference_to_physical_map()", "FEMap");

  // The number of quadrature points.
  const std::size_t n_qp = qp.size();

  // The element type and order to use in
  // the map
  const Order    mapping_order     (elem->default_order());
  const ElemType mapping_elem_type (elem->type());

  // Number of shape functions used to construt the map
  // (Lagrange shape functions are used for mapping)
  const unsigned int n_mapping_shape_functions =
    FE<Dim,LAGRANGE>::n_shape_functions (mapping_elem_type,
                                         mapping_order);

  this->phi_map.resize         (n_mapping_shape_functions);
  if (Dim > 0)
    {
      this->dphidxi_map.resize     (n_mapping_shape_functions);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      this->d2phidxi2_map.resize   (n_mapping_shape_functions);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
    }

  if (Dim > 1)
    {
      this->dphideta_map.resize  (n_mapping_shape_functions);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      this->d2phidxideta_map.resize   (n_mapping_shape_functions);
      this->d2phideta2_map.resize     (n_mapping_shape_functions);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
    }

  if (Dim > 2)
    {
      this->dphidzeta_map.resize (n_mapping_shape_functions);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      this->d2phidxidzeta_map.resize  (n_mapping_shape_functions);
      this->d2phidetadzeta_map.resize (n_mapping_shape_functions);
      this->d2phidzeta2_map.resize    (n_mapping_shape_functions);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
    }


  for (unsigned int i=0; i<n_mapping_shape_functions; i++)
    {
      this->phi_map[i].resize         (n_qp);
      if (Dim > 0)
        {
          this->dphidxi_map[i].resize     (n_qp);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
          this->d2phidxi2_map[i].resize     (n_qp);
          if (Dim > 1)
            {
              this->d2phidxideta_map[i].resize (n_qp);
              this->d2phideta2_map[i].resize (n_qp);
            }
          if (Dim > 2)
            {
              this->d2phidxidzeta_map[i].resize  (n_qp);
              this->d2phidetadzeta_map[i].resize (n_qp);
              this->d2phidzeta2_map[i].resize    (n_qp);
            }
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES

          if (Dim > 1)
            this->dphideta_map[i].resize  (n_qp);

          if (Dim > 2)
            this->dphidzeta_map[i].resize (n_qp);
        }
    }

  // Optimize for the *linear* geometric elements case:
  bool is_linear = elem->is_linear();

  switch (Dim)
    {

      //------------------------------------------------------------
      // 0D
    case 0:
      {
        for (unsigned int i=0; i<n_mapping_shape_functions; i++)
          for (std::size_t p=0; p<n_qp; p++)
            this->phi_map[i][p]      = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i,    qp[p]);

        break;
      }

      //------------------------------------------------------------
      // 1D
    case 1:
      {
        // Compute the value of the mapping shape function i at quadrature point p
        // (Lagrange shape functions are used for mapping)
        if (is_linear)
          {
            for (unsigned int i=0; i<n_mapping_shape_functions; i++)
              {
                this->phi_map[i][0]      = FE<Dim,LAGRANGE>::shape       (mapping_elem_type, mapping_order, i,    qp[0]);
                this->dphidxi_map[i][0]  = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                this->d2phidxi2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                for (std::size_t p=1; p<n_qp; p++)
                  {
                    this->phi_map[i][p]      = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i,    qp[p]);
                    this->dphidxi_map[i][p]  = this->dphidxi_map[i][0];
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                    this->d2phidxi2_map[i][p] = this->d2phidxi2_map[i][0];
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                  }
              }
          }
        else
          for (unsigned int i=0; i<n_mapping_shape_functions; i++)
            for (std::size_t p=0; p<n_qp; p++)
              {
                this->phi_map[i][p]      = FE<Dim,LAGRANGE>::shape       (mapping_elem_type, mapping_order, i,    qp[p]);
                this->dphidxi_map[i][p]  = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                this->d2phidxi2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
              }

        break;
      }
      //------------------------------------------------------------
      // 2D
    case 2:
      {
        // Compute the value of the mapping shape function i at quadrature point p
        // (Lagrange shape functions are used for mapping)
        if (is_linear)
          {
            for (unsigned int i=0; i<n_mapping_shape_functions; i++)
              {
                this->phi_map[i][0]      = FE<Dim,LAGRANGE>::shape       (mapping_elem_type, mapping_order, i,    qp[0]);
                this->dphidxi_map[i][0]  = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]);
                this->dphideta_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[0]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                this->d2phidxi2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]);
                this->d2phidxideta_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 1, qp[0]);
                this->d2phideta2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 2, qp[0]);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                for (std::size_t p=1; p<n_qp; p++)
                  {
                    this->phi_map[i][p]      = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i,    qp[p]);
                    this->dphidxi_map[i][p]  = this->dphidxi_map[i][0];
                    this->dphideta_map[i][p] = this->dphideta_map[i][0];
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                    this->d2phidxi2_map[i][p] = this->d2phidxi2_map[i][0];
                    this->d2phidxideta_map[i][p] = this->d2phidxideta_map[i][0];
                    this->d2phideta2_map[i][p] = this->d2phideta2_map[i][0];
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                  }
              }
          }
        else
          for (unsigned int i=0; i<n_mapping_shape_functions; i++)
            for (std::size_t p=0; p<n_qp; p++)
              {
                this->phi_map[i][p]      = FE<Dim,LAGRANGE>::shape       (mapping_elem_type, mapping_order, i,    qp[p]);
                this->dphidxi_map[i][p]  = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]);
                this->dphideta_map[i][p] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                this->d2phidxi2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]);
                this->d2phidxideta_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]);
                this->d2phideta2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 2, qp[p]);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
              }

        break;
      }

      //------------------------------------------------------------
      // 3D
    case 3:
      {
        // Compute the value of the mapping shape function i at quadrature point p
        // (Lagrange shape functions are used for mapping)
        if (is_linear)
          {
            for (unsigned int i=0; i<n_mapping_shape_functions; i++)
              {
                this->phi_map[i][0]      = FE<Dim,LAGRANGE>::shape       (mapping_elem_type, mapping_order, i,    qp[0]);
                this->dphidxi_map[i][0]  = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]);
                this->dphideta_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[0]);
                this->dphidzeta_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 2, qp[0]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                this->d2phidxi2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]);
                this->d2phidxideta_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 1, qp[0]);
                this->d2phideta2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 2, qp[0]);
                this->d2phidxidzeta_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 3, qp[0]);
                this->d2phidetadzeta_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 4, qp[0]);
                this->d2phidzeta2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 5, qp[0]);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                for (std::size_t p=1; p<n_qp; p++)
                  {
                    this->phi_map[i][p]      = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i,    qp[p]);
                    this->dphidxi_map[i][p]  = this->dphidxi_map[i][0];
                    this->dphideta_map[i][p] = this->dphideta_map[i][0];
                    this->dphidzeta_map[i][p] = this->dphidzeta_map[i][0];
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                    this->d2phidxi2_map[i][p] = this->d2phidxi2_map[i][0];
                    this->d2phidxideta_map[i][p] = this->d2phidxideta_map[i][0];
                    this->d2phideta2_map[i][p] = this->d2phideta2_map[i][0];
                    this->d2phidxidzeta_map[i][p] = this->d2phidxidzeta_map[i][0];
                    this->d2phidetadzeta_map[i][p] = this->d2phidetadzeta_map[i][0];
                    this->d2phidzeta2_map[i][p] = this->d2phidzeta2_map[i][0];
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                  }
              }
          }
        else
          for (unsigned int i=0; i<n_mapping_shape_functions; i++)
            for (std::size_t p=0; p<n_qp; p++)
              {
                this->phi_map[i][p]       = FE<Dim,LAGRANGE>::shape       (mapping_elem_type, mapping_order, i,    qp[p]);
                this->dphidxi_map[i][p]   = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]);
                this->dphideta_map[i][p]  = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]);
                this->dphidzeta_map[i][p] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 2, qp[p]);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
                this->d2phidxi2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]);
                this->d2phidxideta_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]);
                this->d2phideta2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 2, qp[p]);
                this->d2phidxidzeta_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 3, qp[p]);
                this->d2phidetadzeta_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 4, qp[p]);
                this->d2phidzeta2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 5, qp[p]);
#endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
              }

        break;
      }

    default:
      libmesh_error_msg("Invalid Dim = " << Dim);
    }

  // Stop logging the reference->physical map initialization
  STOP_LOG("init_reference_to_physical_map()", "FEMap");
  return;
}
void libMesh::FEMap::print_JxW ( std::ostream &  os) const

Prints the Jacobian times the weight for each quadrature point.

Definition at line 871 of file fe_map.C.

References JxW.

{
  for (unsigned int i=0; i<JxW.size(); ++i)
    os << " [" << i << "]: " <<  JxW[i] << std::endl;
}
void libMesh::FEMap::print_xyz ( std::ostream &  os) const

Prints the spatial location of each quadrature point (on the physical element).

Definition at line 879 of file fe_map.C.

References xyz.

{
  for (unsigned int i=0; i<xyz.size(); ++i)
    os << " [" << i << "]: " << xyz[i];
}
void libMesh::FEMap::resize_quadrature_map_vectors ( const unsigned int  dim,
unsigned int  n_qp 
) [protected]

A utility function for use by compute_*_map

Definition at line 649 of file fe_map.C.

References d2xyzdeta2_map, d2xyzdetadzeta_map, d2xyzdxi2_map, d2xyzdxideta_map, d2xyzdxidzeta_map, d2xyzdzeta2_map, detadx_map, detady_map, detadz_map, dxidx_map, dxidy_map, dxidz_map, dxyzdeta_map, dxyzdxi_map, dxyzdzeta_map, dzetadx_map, dzetady_map, dzetadz_map, jac, JxW, and xyz.

Referenced by compute_affine_map(), compute_map(), and compute_null_map().

{
  // Resize the vectors to hold data at the quadrature points
  xyz.resize(n_qp);
  dxyzdxi_map.resize(n_qp);
  dxidx_map.resize(n_qp);
  dxidy_map.resize(n_qp); // 1D element may live in 2D ...
  dxidz_map.resize(n_qp); // ... or 3D
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
  d2xyzdxi2_map.resize(n_qp);
#endif
  if (dim > 1)
    {
      dxyzdeta_map.resize(n_qp);
      detadx_map.resize(n_qp);
      detady_map.resize(n_qp);
      detadz_map.resize(n_qp);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      d2xyzdxideta_map.resize(n_qp);
      d2xyzdeta2_map.resize(n_qp);
#endif
      if (dim > 2)
        {
          dxyzdzeta_map.resize (n_qp);
          dzetadx_map.resize   (n_qp);
          dzetady_map.resize   (n_qp);
          dzetadz_map.resize   (n_qp);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
          d2xyzdxidzeta_map.resize(n_qp);
          d2xyzdetadzeta_map.resize(n_qp);
          d2xyzdzeta2_map.resize(n_qp);
#endif
        }
    }

  jac.resize(n_qp);
  JxW.resize(n_qp);
}

Member Data Documentation

std::vector<Real> libMesh::FEMap::curvatures [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 735 of file fe_map.h.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and get_curvatures().

std::vector<std::vector<Real> > libMesh::FEMap::d2phideta2_map [protected]

Map for the second derivative, d^2(phi)/d(eta)^2.

Definition at line 668 of file fe_map.h.

Referenced by compute_single_point_map(), get_d2phideta2_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::d2phidetadzeta_map [protected]

Map for the second derivative, d^2(phi)/d(eta)d(zeta).

Definition at line 673 of file fe_map.h.

Referenced by compute_single_point_map(), get_d2phidetadzeta_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::d2phidxi2_map [protected]

Map for the second derivative, d^2(phi)/d(xi)^2.

Definition at line 653 of file fe_map.h.

Referenced by compute_single_point_map(), get_d2phidxi2_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::d2phidxideta_map [protected]

Map for the second derivative, d^2(phi)/d(xi)d(eta).

Definition at line 658 of file fe_map.h.

Referenced by compute_single_point_map(), get_d2phidxideta_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::d2phidxidzeta_map [protected]

Map for the second derivative, d^2(phi)/d(xi)d(zeta).

Definition at line 663 of file fe_map.h.

Referenced by compute_single_point_map(), get_d2phidxidzeta_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::d2phidzeta2_map [protected]

Map for the second derivative, d^2(phi)/d(zeta)^2.

Definition at line 678 of file fe_map.h.

Referenced by compute_single_point_map(), get_d2phidzeta2_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::d2psideta2_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 718 of file fe_map.h.

Referenced by libMesh::FEXYZMap::compute_face_map(), compute_face_map(), get_d2psideta2(), and init_face_shape_functions().

std::vector<std::vector<Real> > libMesh::FEMap::d2psidxi2_map [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 704 of file fe_map.h.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), get_d2psidxi2(), init_edge_shape_functions(), and init_face_shape_functions().

std::vector<std::vector<Real> > libMesh::FEMap::d2psidxideta_map [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 711 of file fe_map.h.

Referenced by libMesh::FEXYZMap::compute_face_map(), compute_face_map(), get_d2psidxideta(), and init_face_shape_functions().

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 562 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_d2xyzdetadzeta(), and resize_quadrature_map_vectors().

std::vector<RealGradient> libMesh::FEMap::d2xyzdxi2_map [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 536 of file fe_map.h.

Referenced by compute_affine_map(), compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), get_d2xyzdxi2(), and resize_quadrature_map_vectors().

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 542 of file fe_map.h.

Referenced by compute_affine_map(), compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), compute_null_map(), compute_single_point_map(), get_d2xyzdxideta(), and resize_quadrature_map_vectors().

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 556 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_d2xyzdxidzeta(), and resize_quadrature_map_vectors().

Vector of second partial derivatives in zeta: d^2(x)/d(zeta)^2

Definition at line 568 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_d2xyzdzeta2(), and resize_quadrature_map_vectors().

std::vector<Real> libMesh::FEMap::detadx_map [protected]

Map for partial derivatives: d(eta)/d(x). Needed for the Jacobian.

Definition at line 595 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_detadx(), and resize_quadrature_map_vectors().

std::vector<Real> libMesh::FEMap::detady_map [protected]

Map for partial derivatives: d(eta)/d(y). Needed for the Jacobian.

Definition at line 601 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_detady(), and resize_quadrature_map_vectors().

std::vector<Real> libMesh::FEMap::detadz_map [protected]

Map for partial derivatives: d(eta)/d(z). Needed for the Jacobian.

Definition at line 607 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_detadz(), and resize_quadrature_map_vectors().

std::vector<std::vector<Real> > libMesh::FEMap::dphideta_map [protected]

Map for the derivative, d(phi)/d(eta).

Definition at line 641 of file fe_map.h.

Referenced by compute_single_point_map(), get_dphideta_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::dphidxi_map [protected]

Map for the derivative, d(phi)/d(xi).

Definition at line 636 of file fe_map.h.

Referenced by compute_single_point_map(), get_dphidxi_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::dphidzeta_map [protected]

Map for the derivative, d(phi)/d(zeta).

Definition at line 646 of file fe_map.h.

Referenced by compute_single_point_map(), get_dphidzeta_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::dpsideta_map [protected]

Map for the derivative of the side function, d(psi)/d(eta).

Definition at line 697 of file fe_map.h.

Referenced by libMesh::FEXYZMap::compute_face_map(), compute_face_map(), get_dpsideta(), and init_face_shape_functions().

std::vector<std::vector<Real> > libMesh::FEMap::dpsidxi_map [protected]

Map for the derivative of the side functions, d(psi)/d(xi).

Definition at line 691 of file fe_map.h.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), get_dpsidxi(), init_edge_shape_functions(), and init_face_shape_functions().

std::vector<Real> libMesh::FEMap::dxidx_map [protected]

Map for partial derivatives: d(xi)/d(x). Needed for the Jacobian.

Definition at line 576 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_dxidx(), and resize_quadrature_map_vectors().

std::vector<Real> libMesh::FEMap::dxidy_map [protected]

Map for partial derivatives: d(xi)/d(y). Needed for the Jacobian.

Definition at line 582 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_dxidy(), and resize_quadrature_map_vectors().

std::vector<Real> libMesh::FEMap::dxidz_map [protected]

Map for partial derivatives: d(xi)/d(z). Needed for the Jacobian.

Definition at line 588 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_dxidz(), and resize_quadrature_map_vectors().

std::vector<RealGradient> libMesh::FEMap::dxyzdzeta_map [protected]

Vector of parital derivatives: d(x)/d(zeta), d(y)/d(zeta), d(z)/d(zeta)

Definition at line 530 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().

std::vector<Real> libMesh::FEMap::dzetadx_map [protected]

Map for partial derivatives: d(zeta)/d(x). Needed for the Jacobian.

Definition at line 614 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_dzetadx(), and resize_quadrature_map_vectors().

std::vector<Real> libMesh::FEMap::dzetady_map [protected]

Map for partial derivatives: d(zeta)/d(y). Needed for the Jacobian.

Definition at line 620 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_dzetady(), and resize_quadrature_map_vectors().

std::vector<Real> libMesh::FEMap::dzetadz_map [protected]

Map for partial derivatives: d(zeta)/d(z). Needed for the Jacobian.

Definition at line 626 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_dzetadz(), and resize_quadrature_map_vectors().

std::vector<Real> libMesh::FEMap::jac [protected]

Jacobian values at quadrature points

Definition at line 740 of file fe_map.h.

Referenced by compute_affine_map(), compute_null_map(), compute_single_point_map(), get_jacobian(), and resize_quadrature_map_vectors().

std::vector<Point> libMesh::FEMap::normals [protected]

Normal vectors on boundary at quadrature points

Definition at line 728 of file fe_map.h.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and get_normals().

std::vector<std::vector<Real> > libMesh::FEMap::phi_map [protected]

Map for the shape function phi.

Definition at line 631 of file fe_map.h.

Referenced by compute_affine_map(), compute_single_point_map(), get_phi_map(), and init_reference_to_physical_map().

std::vector<std::vector<Real> > libMesh::FEMap::psi_map [protected]

Map for the side shape functions, psi.

Definition at line 685 of file fe_map.h.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), get_psi(), init_edge_shape_functions(), and init_face_shape_functions().

std::vector<std::vector<Point> > libMesh::FEMap::tangents [protected]

Tangent vectors on boundary at quadrature points.

Definition at line 723 of file fe_map.h.

Referenced by compute_edge_map(), libMesh::FEXYZMap::compute_face_map(), compute_face_map(), and get_tangents().


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