$extrastylesheet
libMesh::FEAbstract Class Reference

#include <fe_abstract.h>

Inheritance diagram for libMesh::FEAbstract:

List of all members.

Public Member Functions

virtual ~FEAbstract ()
virtual void reinit (const Elem *elem, const std::vector< Point > *const pts=NULL, const std::vector< Real > *const weights=NULL)=0
virtual void reinit (const Elem *elem, const unsigned int side, const Real tolerance=TOLERANCE, const std::vector< Point > *const pts=NULL, const std::vector< Real > *const weights=NULL)=0
virtual void edge_reinit (const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *pts=NULL, const std::vector< Real > *weights=NULL)=0
virtual void side_map (const Elem *elem, const Elem *side, const unsigned int s, const std::vector< Point > &reference_side_points, std::vector< Point > &reference_points)=0
const std::vector< Point > & get_xyz () 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
< Point > > & 
get_tangents () const
const std::vector< Point > & get_normals () const
const std::vector< Real > & get_curvatures () const
virtual void attach_quadrature_rule (QBase *q)=0
virtual unsigned int n_shape_functions () const =0
virtual unsigned int n_quadrature_points () const =0
ElemType get_type () const
unsigned int get_p_level () const
FEType get_fe_type () const
Order get_order () const
virtual FEContinuity get_continuity () const =0
virtual bool is_hierarchic () const =0
FEFamily get_family () const
const FEMapget_fe_map () const
void print_JxW (std::ostream &os) const
virtual void print_phi (std::ostream &os) const =0
virtual void print_dphi (std::ostream &os) const =0
virtual void print_d2phi (std::ostream &os) const =0
void print_xyz (std::ostream &os) const
void print_info (std::ostream &os) const

Static Public Member Functions

static UniquePtr< FEAbstractbuild (const unsigned int dim, const FEType &type)
static bool on_reference_element (const Point &p, const ElemType t, const Real eps=TOLERANCE)
static void get_refspace_nodes (const ElemType t, std::vector< Point > &nodes)
static void compute_node_constraints (NodeConstraints &constraints, const Elem *elem)
static void compute_periodic_node_constraints (NodeConstraints &constraints, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const Elem *elem)
static std::string get_info ()
static void print_info (std::ostream &out=libMesh::out)
static unsigned int n_objects ()
static void enable_print_counter_info ()
static void disable_print_counter_info ()

Protected Types

typedef std::map< std::string,
std::pair< unsigned int,
unsigned int > > 
Counts

Protected Member Functions

 FEAbstract (const unsigned int dim, const FEType &fet)
virtual void compute_shape_functions (const Elem *, const std::vector< Point > &)=0
virtual bool shapes_need_reinit () const =0
void increment_constructor_count (const std::string &name)
void increment_destructor_count (const std::string &name)

Protected Attributes

UniquePtr< FEMap_fe_map
const unsigned int dim
bool calculations_started
bool calculate_phi
bool calculate_dphi
bool calculate_d2phi
bool calculate_curl_phi
bool calculate_div_phi
bool calculate_dphiref
const FEType fe_type
ElemType elem_type
unsigned int _p_level
QBaseqrule
bool shapes_on_quadrature

Static Protected Attributes

static Counts _counts
static Threads::atomic
< unsigned int > 
_n_objects
static Threads::spin_mutex _mutex
static bool _enable_print_counter = true

Friends

std::ostream & operator<< (std::ostream &os, const FEAbstract &fe)

Detailed Description

This class forms the foundation from which generic finite elements may be derived. In the current implementation the templated derived class FE offers a wide variety of commonly used finite element concepts. Check there for details. Use the FEAbstract::build() method to create an object of any of the derived classes. Note that the amount of virtual members is kept to a minimum, and the sophisticated template scheme of FE is quite likely to offer acceptably fast code.

All calls to static members of the FE classes should be requested through the FEInterface. This interface class offers sort-of runtime polymorphism for the templated finite element classes. Even internal library classes, like DofMap, request the number of dof's through this interface class. Note that this also enables the co-existence of various element-based schemes. This class is well 'at the heart' of the library, so things in here should better remain unchanged.

Author:
Benjamin S. Kirk, 2002

Definition at line 97 of file fe_abstract.h.


Member Typedef Documentation

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited]

Data structure to log the information. The log is identified by the class name.

Definition at line 113 of file reference_counter.h.


Constructor & Destructor Documentation

libMesh::FEAbstract::FEAbstract ( const unsigned int  dim,
const FEType fet 
) [inline, protected]

Constructor. Optionally initializes required data structures. Protected so that this base class cannot be explicitly instantiated.

Definition at line 599 of file fe_abstract.h.

libMesh::FEAbstract::~FEAbstract ( ) [inline, virtual]

Destructor.

Definition at line 620 of file fe_abstract.h.

{
}

Member Function Documentation

UniquePtr< FEAbstract > libMesh::FEAbstract::build ( const unsigned int  dim,
const FEType type 
) [static]

Builds a specific finite element type. A UniquePtr<FEAbstract> is returned to prevent a memory leak. This way the user need not remember to delete the object.

Reimplemented in libMesh::FEGenericBase< OutputType >, libMesh::FEGenericBase< OutputType >, libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.

Definition at line 44 of file fe_abstract.C.

References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::L2_HIERARCHIC, libMesh::L2_LAGRANGE, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::NEDELEC_ONE, libMesh::SCALAR, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.

Referenced by libMesh::DGFEMContext::DGFEMContext(), libMesh::FEMContext::FEMContext(), and libMesh::DofMap::use_coupled_neighbor_dofs().

{
  switch (dim)
    {
      // 0D
    case 0:
      {
        switch (fet.family)
          {
          case CLOUGH:
            return UniquePtr<FEAbstract>(new FE<0,CLOUGH>(fet));

          case HERMITE:
            return UniquePtr<FEAbstract>(new FE<0,HERMITE>(fet));

          case LAGRANGE:
            return UniquePtr<FEAbstract>(new FE<0,LAGRANGE>(fet));

          case LAGRANGE_VEC:
            return UniquePtr<FEAbstract>(new FE<0,LAGRANGE_VEC>(fet));

          case L2_LAGRANGE:
            return UniquePtr<FEAbstract>(new FE<0,L2_LAGRANGE>(fet));

          case HIERARCHIC:
            return UniquePtr<FEAbstract>(new FE<0,HIERARCHIC>(fet));

          case L2_HIERARCHIC:
            return UniquePtr<FEAbstract>(new FE<0,L2_HIERARCHIC>(fet));

          case MONOMIAL:
            return UniquePtr<FEAbstract>(new FE<0,MONOMIAL>(fet));

#ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
          case SZABAB:
            return UniquePtr<FEAbstract>(new FE<0,SZABAB>(fet));

          case BERNSTEIN:
            return UniquePtr<FEAbstract>(new FE<0,BERNSTEIN>(fet));
#endif

          case XYZ:
            return UniquePtr<FEAbstract>(new FEXYZ<0>(fet));

          case SCALAR:
            return UniquePtr<FEAbstract>(new FEScalar<0>(fet));

          default:
            libmesh_error_msg("ERROR: Bad FEType.family= " << fet.family);
          }
      }
      // 1D
    case 1:
      {
        switch (fet.family)
          {
          case CLOUGH:
            return UniquePtr<FEAbstract>(new FE<1,CLOUGH>(fet));

          case HERMITE:
            return UniquePtr<FEAbstract>(new FE<1,HERMITE>(fet));

          case LAGRANGE:
            return UniquePtr<FEAbstract>(new FE<1,LAGRANGE>(fet));

          case LAGRANGE_VEC:
            return UniquePtr<FEAbstract>(new FE<1,LAGRANGE_VEC>(fet));

          case L2_LAGRANGE:
            return UniquePtr<FEAbstract>(new FE<1,L2_LAGRANGE>(fet));

          case HIERARCHIC:
            return UniquePtr<FEAbstract>(new FE<1,HIERARCHIC>(fet));

          case L2_HIERARCHIC:
            return UniquePtr<FEAbstract>(new FE<1,L2_HIERARCHIC>(fet));

          case MONOMIAL:
            return UniquePtr<FEAbstract>(new FE<1,MONOMIAL>(fet));

#ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
          case SZABAB:
            return UniquePtr<FEAbstract>(new FE<1,SZABAB>(fet));

          case BERNSTEIN:
            return UniquePtr<FEAbstract>(new FE<1,BERNSTEIN>(fet));
#endif

          case XYZ:
            return UniquePtr<FEAbstract>(new FEXYZ<1>(fet));

          case SCALAR:
            return UniquePtr<FEAbstract>(new FEScalar<1>(fet));

          default:
            libmesh_error_msg("ERROR: Bad FEType.family= " << fet.family);
          }
      }


      // 2D
    case 2:
      {
        switch (fet.family)
          {
          case CLOUGH:
            return UniquePtr<FEAbstract>(new FE<2,CLOUGH>(fet));

          case HERMITE:
            return UniquePtr<FEAbstract>(new FE<2,HERMITE>(fet));

          case LAGRANGE:
            return UniquePtr<FEAbstract>(new FE<2,LAGRANGE>(fet));

          case LAGRANGE_VEC:
            return UniquePtr<FEAbstract>(new FE<2,LAGRANGE_VEC>(fet));

          case L2_LAGRANGE:
            return UniquePtr<FEAbstract>(new FE<2,L2_LAGRANGE>(fet));

          case HIERARCHIC:
            return UniquePtr<FEAbstract>(new FE<2,HIERARCHIC>(fet));

          case L2_HIERARCHIC:
            return UniquePtr<FEAbstract>(new FE<2,L2_HIERARCHIC>(fet));

          case MONOMIAL:
            return UniquePtr<FEAbstract>(new FE<2,MONOMIAL>(fet));

#ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
          case SZABAB:
            return UniquePtr<FEAbstract>(new FE<2,SZABAB>(fet));

          case BERNSTEIN:
            return UniquePtr<FEAbstract>(new FE<2,BERNSTEIN>(fet));
#endif

          case XYZ:
            return UniquePtr<FEAbstract>(new FEXYZ<2>(fet));

          case SCALAR:
            return UniquePtr<FEAbstract>(new FEScalar<2>(fet));

          case NEDELEC_ONE:
            return UniquePtr<FEAbstract>(new FENedelecOne<2>(fet));

          case SUBDIVISION:
            return UniquePtr<FEAbstract>(new FESubdivision(fet));

          default:
            libmesh_error_msg("ERROR: Bad FEType.family= " << fet.family);
          }
      }


      // 3D
    case 3:
      {
        switch (fet.family)
          {
          case CLOUGH:
            libmesh_error_msg("ERROR: Clough-Tocher elements currently only support 1D and 2D");

          case HERMITE:
            return UniquePtr<FEAbstract>(new FE<3,HERMITE>(fet));

          case LAGRANGE:
            return UniquePtr<FEAbstract>(new FE<3,LAGRANGE>(fet));

          case LAGRANGE_VEC:
            return UniquePtr<FEAbstract>(new FE<3,LAGRANGE_VEC>(fet));

          case L2_LAGRANGE:
            return UniquePtr<FEAbstract>(new FE<3,L2_LAGRANGE>(fet));

          case HIERARCHIC:
            return UniquePtr<FEAbstract>(new FE<3,HIERARCHIC>(fet));

          case L2_HIERARCHIC:
            return UniquePtr<FEAbstract>(new FE<3,L2_HIERARCHIC>(fet));

          case MONOMIAL:
            return UniquePtr<FEAbstract>(new FE<3,MONOMIAL>(fet));

#ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
          case SZABAB:
            return UniquePtr<FEAbstract>(new FE<3,SZABAB>(fet));

          case BERNSTEIN:
            return UniquePtr<FEAbstract>(new FE<3,BERNSTEIN>(fet));
#endif

          case XYZ:
            return UniquePtr<FEAbstract>(new FEXYZ<3>(fet));

          case SCALAR:
            return UniquePtr<FEAbstract>(new FEScalar<3>(fet));

          case NEDELEC_ONE:
            return UniquePtr<FEAbstract>(new FENedelecOne<3>(fet));

          default:
            libmesh_error_msg("ERROR: Bad FEType.family= " << fet.family);
          }
      }

    default:
      libmesh_error_msg("Invalid dimension dim = " << dim);
    }

  libmesh_error_msg("We'll never get here!");
  return UniquePtr<FEAbstract>();
}
void libMesh::FEAbstract::compute_node_constraints ( NodeConstraints constraints,
const Elem elem 
) [static]

Computes the nodal constraint contributions (for non-conforming adapted meshes), using Lagrange geometry

Definition at line 788 of file fe_abstract.C.

References std::abs(), libMesh::Elem::build_side(), libMesh::Elem::default_order(), libMesh::Elem::dim(), fe_type, libMesh::FEInterface::inverse_map(), libMesh::LAGRANGE, libMesh::Elem::level(), libMesh::libmesh_assert(), libMesh::FEInterface::n_dofs(), libMesh::Elem::n_sides(), libMesh::Elem::neighbor(), libMesh::Elem::parent(), libMesh::Real, libMesh::remote_elem, libMesh::FEInterface::shape(), libMesh::Threads::spin_mtx, and libMesh::Elem::subactive().

{
  libmesh_assert(elem);

  const unsigned int Dim = elem->dim();

  // Only constrain elements in 2,3D.
  if (Dim == 1)
    return;

  // Only constrain active and ancestor elements
  if (elem->subactive())
    return;

  // We currently always use LAGRANGE mappings for geometry
  const FEType fe_type(elem->default_order(), LAGRANGE);

  std::vector<const Node*> my_nodes, parent_nodes;

  // Look at the element faces.  Check to see if we need to
  // build constraints.
  for (unsigned int s=0; s<elem->n_sides(); s++)
    if (elem->neighbor(s) != NULL &&
        elem->neighbor(s) != remote_elem)
      if (elem->neighbor(s)->level() < elem->level()) // constrain dofs shared between
        {                                                     // this element and ones coarser
          // than this element.
          // Get pointers to the elements of interest and its parent.
          const Elem* parent = elem->parent();

          // This can't happen...  Only level-0 elements have NULL
          // parents, and no level-0 elements can be at a higher
          // level than their neighbors!
          libmesh_assert(parent);

          const UniquePtr<Elem> my_side     (elem->build_side(s));
          const UniquePtr<Elem> parent_side (parent->build_side(s));

          const unsigned int n_side_nodes = my_side->n_nodes();

          my_nodes.clear();
          my_nodes.reserve (n_side_nodes);
          parent_nodes.clear();
          parent_nodes.reserve (n_side_nodes);

          for (unsigned int n=0; n != n_side_nodes; ++n)
            my_nodes.push_back(my_side->get_node(n));

          for (unsigned int n=0; n != n_side_nodes; ++n)
            parent_nodes.push_back(parent_side->get_node(n));

          for (unsigned int my_side_n=0;
               my_side_n < n_side_nodes;
               my_side_n++)
            {
              libmesh_assert_less (my_side_n, FEInterface::n_dofs(Dim-1, fe_type, my_side->type()));

              const Node* my_node = my_nodes[my_side_n];

              // The support point of the DOF
              const Point& support_point = *my_node;

              // Figure out where my node lies on their reference element.
              const Point mapped_point = FEInterface::inverse_map(Dim-1, fe_type,
                                                                  parent_side.get(),
                                                                  support_point);

              // Compute the parent's side shape function values.
              for (unsigned int their_side_n=0;
                   their_side_n < n_side_nodes;
                   their_side_n++)
                {
                  libmesh_assert_less (their_side_n, FEInterface::n_dofs(Dim-1, fe_type, parent_side->type()));

                  const Node* their_node = parent_nodes[their_side_n];
                  libmesh_assert(their_node);

                  const Real their_value = FEInterface::shape(Dim-1,
                                                              fe_type,
                                                              parent_side->type(),
                                                              their_side_n,
                                                              mapped_point);

                  const Real their_mag = std::abs(their_value);
#ifdef DEBUG
                  // Protect for the case u_i ~= u_j,
                  // in which case i better equal j.
                  if (their_mag > 0.999)
                    {
                      libmesh_assert_equal_to (my_node, their_node);
                      libmesh_assert_less (std::abs(their_value - 1.), 0.001);
                    }
                  else
#endif
                    // To make nodal constraints useful for constructing
                    // sparsity patterns faster, we need to get EVERY
                    // POSSIBLE constraint coupling identified, even if
                    // there is no coupling in the isoparametric
                    // Lagrange case.
                    if (their_mag < 1.e-5)
                      {
                        // since we may be running this method concurrently
                        // on multiple threads we need to acquire a lock
                        // before modifying the shared constraint_row object.
                        Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);

                        // A reference to the constraint row.
                        NodeConstraintRow& constraint_row = constraints[my_node].first;

                        constraint_row.insert(std::make_pair (their_node,
                                                              0.));
                      }
                  // To get nodal coordinate constraints right, only
                  // add non-zero and non-identity values for Lagrange
                  // basis functions.
                    else // (1.e-5 <= their_mag <= .999)
                      {
                        // since we may be running this method concurrently
                        // on multiple threads we need to acquire a lock
                        // before modifying the shared constraint_row object.
                        Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);

                        // A reference to the constraint row.
                        NodeConstraintRow& constraint_row = constraints[my_node].first;

                        constraint_row.insert(std::make_pair (their_node,
                                                              their_value));
                      }
                }
            }
        }
}
void libMesh::FEAbstract::compute_periodic_node_constraints ( NodeConstraints constraints,
const PeriodicBoundaries boundaries,
const MeshBase mesh,
const PointLocatorBase point_locator,
const Elem elem 
) [static]

Computes the node position constraint equation contributions (for meshes with periodic boundary conditions)

Definition at line 931 of file fe_abstract.C.

References libMesh::Elem::active(), libMesh::PeriodicBoundaries::boundary(), libMesh::BoundaryInfo::boundary_ids(), libMesh::Elem::build_side(), libMesh::Elem::default_order(), libMesh::Elem::dim(), fe_type, libMesh::MeshBase::get_boundary_info(), libMesh::PeriodicBoundaryBase::get_corresponding_pos(), libMesh::invalid_uint, libMesh::FEInterface::inverse_map(), libMesh::LAGRANGE, libMesh::Elem::level(), libMesh::libmesh_assert(), libMesh::FEInterface::n_dofs(), libMesh::Elem::n_sides(), libMesh::PeriodicBoundaries::neighbor(), libMesh::Elem::neighbor(), libMesh::PeriodicBoundaryBase::pairedboundary, libMesh::Real, libMesh::FEInterface::shape(), libMesh::BoundaryInfo::side_with_boundary_id(), and libMesh::Threads::spin_mtx.

{
  // Only bother if we truly have periodic boundaries
  if (boundaries.empty())
    return;

  libmesh_assert(elem);

  // Only constrain active elements with this method
  if (!elem->active())
    return;

  const unsigned int Dim = elem->dim();

  // We currently always use LAGRANGE mappings for geometry
  const FEType fe_type(elem->default_order(), LAGRANGE);

  std::vector<const Node*> my_nodes, neigh_nodes;

  // Look at the element faces.  Check to see if we need to
  // build constraints.
  for (unsigned short int s=0; s<elem->n_sides(); s++)
    {
      if (elem->neighbor(s))
        continue;

      const std::vector<boundary_id_type>& bc_ids =
        mesh.get_boundary_info().boundary_ids (elem, s);
      for (std::vector<boundary_id_type>::const_iterator id_it=bc_ids.begin(); id_it!=bc_ids.end(); ++id_it)
        {
          const boundary_id_type boundary_id = *id_it;
          const PeriodicBoundaryBase *periodic = boundaries.boundary(boundary_id);
          if (periodic)
            {
              libmesh_assert(point_locator);

              // Get pointers to the element's neighbor.
              const Elem* neigh = boundaries.neighbor(boundary_id, *point_locator, elem, s);

              // h refinement constraints:
              // constrain dofs shared between
              // this element and ones as coarse
              // as or coarser than this element.
              if (neigh->level() <= elem->level())
                {
                  unsigned int s_neigh =
                    mesh.get_boundary_info().side_with_boundary_id(neigh, periodic->pairedboundary);
                  libmesh_assert_not_equal_to (s_neigh, libMesh::invalid_uint);

#ifdef LIBMESH_ENABLE_AMR
                  libmesh_assert(neigh->active());
#endif // #ifdef LIBMESH_ENABLE_AMR

                  const UniquePtr<Elem> my_side    (elem->build_side(s));
                  const UniquePtr<Elem> neigh_side (neigh->build_side(s_neigh));

                  const unsigned int n_side_nodes = my_side->n_nodes();

                  my_nodes.clear();
                  my_nodes.reserve (n_side_nodes);
                  neigh_nodes.clear();
                  neigh_nodes.reserve (n_side_nodes);

                  for (unsigned int n=0; n != n_side_nodes; ++n)
                    my_nodes.push_back(my_side->get_node(n));

                  for (unsigned int n=0; n != n_side_nodes; ++n)
                    neigh_nodes.push_back(neigh_side->get_node(n));

                  // Make sure we're not adding recursive constraints
                  // due to the redundancy in the way we add periodic
                  // boundary constraints, or adding constraints to
                  // nodes that already have AMR constraints
                  std::vector<bool> skip_constraint(n_side_nodes, false);

                  for (unsigned int my_side_n=0;
                       my_side_n < n_side_nodes;
                       my_side_n++)
                    {
                      libmesh_assert_less (my_side_n, FEInterface::n_dofs(Dim-1, fe_type, my_side->type()));

                      const Node* my_node = my_nodes[my_side_n];

                      // Figure out where my node lies on their reference element.
                      const Point neigh_point = periodic->get_corresponding_pos(*my_node);

                      const Point mapped_point = FEInterface::inverse_map(Dim-1, fe_type,
                                                                          neigh_side.get(),
                                                                          neigh_point);

                      // If we've already got a constraint on this
                      // node, then the periodic constraint is
                      // redundant
                      {
                        Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);

                        if (constraints.count(my_node))
                          {
                            skip_constraint[my_side_n] = true;
                            continue;
                          }
                      }

                      // Compute the neighbors's side shape function values.
                      for (unsigned int their_side_n=0;
                           their_side_n < n_side_nodes;
                           their_side_n++)
                        {
                          libmesh_assert_less (their_side_n, FEInterface::n_dofs(Dim-1, fe_type, neigh_side->type()));

                          const Node* their_node = neigh_nodes[their_side_n];

                          // If there's a constraint on an opposing node,
                          // we need to see if it's constrained by
                          // *our side* making any periodic constraint
                          // on us recursive
                          {
                            Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);

                            if (!constraints.count(their_node))
                              continue;

                            const NodeConstraintRow& their_constraint_row =
                              constraints[their_node].first;

                            for (unsigned int orig_side_n=0;
                                 orig_side_n < n_side_nodes;
                                 orig_side_n++)
                              {
                                libmesh_assert_less (orig_side_n, FEInterface::n_dofs(Dim-1, fe_type, my_side->type()));

                                const Node* orig_node = my_nodes[orig_side_n];

                                if (their_constraint_row.count(orig_node))
                                  skip_constraint[orig_side_n] = true;
                              }
                          }
                        }
                    }
                  for (unsigned int my_side_n=0;
                       my_side_n < n_side_nodes;
                       my_side_n++)
                    {
                      libmesh_assert_less (my_side_n, FEInterface::n_dofs(Dim-1, fe_type, my_side->type()));

                      if (skip_constraint[my_side_n])
                        continue;

                      const Node* my_node = my_nodes[my_side_n];

                      // Figure out where my node lies on their reference element.
                      const Point neigh_point = periodic->get_corresponding_pos(*my_node);

                      // Figure out where my node lies on their reference element.
                      const Point mapped_point = FEInterface::inverse_map(Dim-1, fe_type,
                                                                          neigh_side.get(),
                                                                          neigh_point);

                      for (unsigned int their_side_n=0;
                           their_side_n < n_side_nodes;
                           their_side_n++)
                        {
                          libmesh_assert_less (their_side_n, FEInterface::n_dofs(Dim-1, fe_type, neigh_side->type()));

                          const Node* their_node = neigh_nodes[their_side_n];
                          libmesh_assert(their_node);

                          const Real their_value = FEInterface::shape(Dim-1,
                                                                      fe_type,
                                                                      neigh_side->type(),
                                                                      their_side_n,
                                                                      mapped_point);

                          // since we may be running this method concurrently
                          // on multiple threads we need to acquire a lock
                          // before modifying the shared constraint_row object.
                          {
                            Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);

                            NodeConstraintRow& constraint_row =
                              constraints[my_node].first;

                            constraint_row.insert(std::make_pair(their_node,
                                                                 their_value));
                          }
                        }
                    }
                }
            }
        }
    }
}
virtual void libMesh::FEAbstract::compute_shape_functions ( const Elem ,
const std::vector< Point > &   
) [protected, pure virtual]

After having updated the jacobian and the transformation from local to global coordinates in FEMap::compute_map(), the first derivatives of the shape functions are transformed to global coordinates, giving dphi, dphidx, dphidy, and dphidz. This method should rarely be re-defined in derived classes, but still should be usable for children. Therefore, keep it protected. This needs to be implemented in the derived class since this function depends on whether the shape functions are vector-valued or not.

Implemented in libMesh::FEXYZ< Dim >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.

virtual void libMesh::FEAbstract::edge_reinit ( const Elem elem,
const unsigned int  edge,
const Real  tolerance = TOLERANCE,
const std::vector< Point > *  pts = NULL,
const std::vector< Real > *  weights = NULL 
) [pure virtual]

Reinitializes all the physical element-dependent data based on the edge of the element elem. The tolerance paremeter is passed to the involved call to inverse_map(). By default the element data are computed at the quadrature points specified by the quadrature rule qrule, but any set of points on the reference edge element may be specified in the optional argument pts.

Implemented in libMesh::FE< Dim, T >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, and libMesh::FE< Dim, LAGRANGE_VEC >.

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

{
  _enable_print_counter = true;
  return;
}
virtual FEContinuity libMesh::FEAbstract::get_continuity ( ) const [pure virtual]
Returns:
the continuity level of the finite element.

Implemented in libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, libMesh::FE< Dim, LAGRANGE_VEC >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, and libMesh::FE< Dim, T >.

Referenced by libMesh::ProjectFEMSolution::operator()().

const std::vector<Real>& libMesh::FEAbstract::get_curvatures ( ) const [inline]
Returns:
the curvatures for use in face integration.

Definition at line 380 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_curvatures();}
const std::vector<RealGradient>& libMesh::FEAbstract::get_d2xyzdeta2 ( ) const [inline]
Returns:
the second partial derivatives in eta.

Definition at line 267 of file fe_abstract.h.

References _fe_map.

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

Definition at line 297 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_d2xyzdetadzeta(); }
const std::vector<RealGradient>& libMesh::FEAbstract::get_d2xyzdxi2 ( ) const [inline]
Returns:
the second partial derivatives in xi.

Definition at line 261 of file fe_abstract.h.

References _fe_map.

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

Definition at line 283 of file fe_abstract.h.

References _fe_map.

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

Definition at line 291 of file fe_abstract.h.

References _fe_map.

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

Definition at line 275 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_d2xyzdzeta2(); }
const std::vector<Real>& libMesh::FEAbstract::get_detadx ( ) const [inline]
Returns:
the deta/dx entry in the transformation matrix from physical to local coordinates.

Definition at line 327 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_detadx(); }
const std::vector<Real>& libMesh::FEAbstract::get_detady ( ) const [inline]
Returns:
the deta/dy entry in the transformation matrix from physical to local coordinates.

Definition at line 334 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_detady(); }
const std::vector<Real>& libMesh::FEAbstract::get_detadz ( ) const [inline]
Returns:
the deta/dz entry in the transformation matrix from physical to local coordinates.

Definition at line 341 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_detadz(); }
const std::vector<Real>& libMesh::FEAbstract::get_dxidx ( ) const [inline]
Returns:
the dxi/dx entry in the transformation matrix from physical to local coordinates.

Definition at line 306 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_dxidx(); }
const std::vector<Real>& libMesh::FEAbstract::get_dxidy ( ) const [inline]
Returns:
the dxi/dy entry in the transformation matrix from physical to local coordinates.

Definition at line 313 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_dxidy(); }
const std::vector<Real>& libMesh::FEAbstract::get_dxidz ( ) const [inline]
Returns:
the dxi/dz entry in the transformation matrix from physical to local coordinates.

Definition at line 320 of file fe_abstract.h.

References _fe_map.

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

Definition at line 248 of file fe_abstract.h.

References _fe_map.

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

Definition at line 241 of file fe_abstract.h.

References _fe_map.

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

Definition at line 255 of file fe_abstract.h.

References _fe_map.

  { return _fe_map->get_dxyzdzeta(); }
const std::vector<Real>& libMesh::FEAbstract::get_dzetadx ( ) const [inline]
Returns:
the dzeta/dx entry in the transformation matrix from physical to local coordinates.

Definition at line 348 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_dzetadx(); }
const std::vector<Real>& libMesh::FEAbstract::get_dzetady ( ) const [inline]
Returns:
the dzeta/dy entry in the transformation matrix from physical to local coordinates.

Definition at line 355 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_dzetady(); }
const std::vector<Real>& libMesh::FEAbstract::get_dzetadz ( ) const [inline]
Returns:
the dzeta/dz entry in the transformation matrix from physical to local coordinates.

Definition at line 362 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_dzetadz(); }
Returns:
the finite element family of this element.

Definition at line 439 of file fe_abstract.h.

References libMesh::FEType::family, and fe_type.

Referenced by libMesh::FE< Dim, T >::FE().

{ return fe_type.family; }
Returns:
the FE Type (approximation order and family) of the finite element.

Definition at line 418 of file fe_abstract.h.

References fe_type.

Referenced by libMesh::FEMContext::build_new_fe(), libMesh::H1FETransformation< OutputShape >::map_phi(), libMesh::HCurlFETransformation< OutputShape >::map_phi(), and libMesh::ProjectFEMSolution::operator()().

{ return fe_type; }
std::string libMesh::ReferenceCounter::get_info ( ) [static, inherited]

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)

  std::ostringstream oss;

  oss << '\n'
      << " ---------------------------------------------------------------------------- \n"
      << "| Reference count information                                                |\n"
      << " ---------------------------------------------------------------------------- \n";

  for (Counts::iterator it = _counts.begin();
       it != _counts.end(); ++it)
    {
      const std::string name(it->first);
      const unsigned int creations    = it->second.first;
      const unsigned int destructions = it->second.second;

      oss << "| " << name << " reference count information:\n"
          << "|  Creations:    " << creations    << '\n'
          << "|  Destructions: " << destructions << '\n';
    }

  oss << " ---------------------------------------------------------------------------- \n";

  return oss.str();

#else

  return "";

#endif
}
const std::vector<Real>& libMesh::FEAbstract::get_JxW ( ) const [inline]
Returns:
the element Jacobian times the quadrature weight for each quadrature point.

Definition at line 234 of file fe_abstract.h.

References _fe_map.

Referenced by libMesh::ExactSolution::_compute_error(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::FEMSystem::init_context(), and libMesh::ProjectFEMSolution::operator()().

  { return this->_fe_map->get_JxW(); }
const std::vector<Point>& libMesh::FEAbstract::get_normals ( ) const [inline]
Returns:
the normal vectors for face integration.

Definition at line 374 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_normals(); }
Returns:
the approximation order of the finite element.

Definition at line 423 of file fe_abstract.h.

References _p_level, fe_type, and libMesh::FEType::order.

{ return static_cast<Order>(fe_type.order + _p_level); }
unsigned int libMesh::FEAbstract::get_p_level ( ) const [inline]
Returns:
the p refinement level that the current shape functions have been calculated for.

Definition at line 413 of file fe_abstract.h.

References _p_level.

{ return _p_level; }
void libMesh::FEAbstract::get_refspace_nodes ( const ElemType  t,
std::vector< Point > &  nodes 
) [static]

returns the reference space nodes coordinates given the element type

Definition at line 259 of file fe_abstract.C.

References libMesh::EDGE2, libMesh::EDGE3, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::TET10, libMesh::TET4, libMesh::TRI3, and libMesh::TRI6.

{
  switch(itemType)
    {
    case EDGE2:
      {
        nodes.resize(2);
        nodes[0] = Point (-1.,0.,0.);
        nodes[1] = Point (1.,0.,0.);
        return;
      }
    case EDGE3:
      {
        nodes.resize(3);
        nodes[0] = Point (-1.,0.,0.);
        nodes[1] = Point (1.,0.,0.);
        nodes[2] = Point (0.,0.,0.);
        return;
      }
    case TRI3:
      {
        nodes.resize(3);
        nodes[0] = Point (0.,0.,0.);
        nodes[1] = Point (1.,0.,0.);
        nodes[2] = Point (0.,1.,0.);
        return;
      }
    case TRI6:
      {
        nodes.resize(6);
        nodes[0] = Point (0.,0.,0.);
        nodes[1] = Point (1.,0.,0.);
        nodes[2] = Point (0.,1.,0.);
        nodes[3] = Point (.5,0.,0.);
        nodes[4] = Point (.5,.5,0.);
        nodes[5] = Point (0.,.5,0.);
        return;
      }
    case QUAD4:
      {
        nodes.resize(4);
        nodes[0] = Point (-1.,-1.,0.);
        nodes[1] = Point (1.,-1.,0.);
        nodes[2] = Point (1.,1.,0.);
        nodes[3] = Point (-1.,1.,0.);
        return;
      }
    case QUAD8:
      {
        nodes.resize(8);
        nodes[0] = Point (-1.,-1.,0.);
        nodes[1] = Point (1.,-1.,0.);
        nodes[2] = Point (1.,1.,0.);
        nodes[3] = Point (-1.,1.,0.);
        nodes[4] = Point (0.,-1.,0.);
        nodes[5] = Point (1.,0.,0.);
        nodes[6] = Point (0.,1.,0.);
        nodes[7] = Point (-1.,0.,0.);
        return;
      }
    case QUAD9:
      {
        nodes.resize(9);
        nodes[0] = Point (-1.,-1.,0.);
        nodes[1] = Point (1.,-1.,0.);
        nodes[2] = Point (1.,1.,0.);
        nodes[3] = Point (-1.,1.,0.);
        nodes[4] = Point (0.,-1.,0.);
        nodes[5] = Point (1.,0.,0.);
        nodes[6] = Point (0.,1.,0.);
        nodes[7] = Point (-1.,0.,0.);
        nodes[8] = Point (0.,0.,0.);
        return;
      }
    case TET4:
      {
        nodes.resize(4);
        nodes[0] = Point (0.,0.,0.);
        nodes[1] = Point (1.,0.,0.);
        nodes[2] = Point (0.,1.,0.);
        nodes[3] = Point (0.,0.,1.);
        return;
      }
    case TET10:
      {
        nodes.resize(10);
        nodes[0] = Point (0.,0.,0.);
        nodes[1] = Point (1.,0.,0.);
        nodes[2] = Point (0.,1.,0.);
        nodes[3] = Point (0.,0.,1.);
        nodes[4] = Point (.5,0.,0.);
        nodes[5] = Point (.5,.5,0.);
        nodes[6] = Point (0.,.5,0.);
        nodes[7] = Point (0.,0.,.5);
        nodes[8] = Point (.5,0.,.5);
        nodes[9] = Point (0.,.5,.5);
        return;
      }
    case HEX8:
      {
        nodes.resize(8);
        nodes[0] = Point (-1.,-1.,-1.);
        nodes[1] = Point (1.,-1.,-1.);
        nodes[2] = Point (1.,1.,-1.);
        nodes[3] = Point (-1.,1.,-1.);
        nodes[4] = Point (-1.,-1.,1.);
        nodes[5] = Point (1.,-1.,1.);
        nodes[6] = Point (1.,1.,1.);
        nodes[7] = Point (-1.,1.,1.);
        return;
      }
    case HEX20:
      {
        nodes.resize(20);
        nodes[0] = Point (-1.,-1.,-1.);
        nodes[1] = Point (1.,-1.,-1.);
        nodes[2] = Point (1.,1.,-1.);
        nodes[3] = Point (-1.,1.,-1.);
        nodes[4] = Point (-1.,-1.,1.);
        nodes[5] = Point (1.,-1.,1.);
        nodes[6] = Point (1.,1.,1.);
        nodes[7] = Point (-1.,1.,1.);
        nodes[8] = Point (0.,-1.,-1.);
        nodes[9] = Point (1.,0.,-1.);
        nodes[10] = Point (0.,1.,-1.);
        nodes[11] = Point (-1.,0.,-1.);
        nodes[12] = Point (-1.,-1.,0.);
        nodes[13] = Point (1.,-1.,0.);
        nodes[14] = Point (1.,1.,0.);
        nodes[15] = Point (-1.,1.,0.);
        nodes[16] = Point (0.,-1.,1.);
        nodes[17] = Point (1.,0.,1.);
        nodes[18] = Point (0.,1.,1.);
        nodes[19] = Point (-1.,0.,1.);
        return;
      }
    case HEX27:
      {
        nodes.resize(27);
        nodes[0] = Point (-1.,-1.,-1.);
        nodes[1] = Point (1.,-1.,-1.);
        nodes[2] = Point (1.,1.,-1.);
        nodes[3] = Point (-1.,1.,-1.);
        nodes[4] = Point (-1.,-1.,1.);
        nodes[5] = Point (1.,-1.,1.);
        nodes[6] = Point (1.,1.,1.);
        nodes[7] = Point (-1.,1.,1.);
        nodes[8] = Point (0.,-1.,-1.);
        nodes[9] = Point (1.,0.,-1.);
        nodes[10] = Point (0.,1.,-1.);
        nodes[11] = Point (-1.,0.,-1.);
        nodes[12] = Point (-1.,-1.,0.);
        nodes[13] = Point (1.,-1.,0.);
        nodes[14] = Point (1.,1.,0.);
        nodes[15] = Point (-1.,1.,0.);
        nodes[16] = Point (0.,-1.,1.);
        nodes[17] = Point (1.,0.,1.);
        nodes[18] = Point (0.,1.,1.);
        nodes[19] = Point (-1.,0.,1.);
        nodes[20] = Point (0.,0.,-1.);
        nodes[21] = Point (0.,-1.,0.);
        nodes[22] = Point (1.,0.,0.);
        nodes[23] = Point (0.,1.,0.);
        nodes[24] = Point (-1.,0.,0.);
        nodes[25] = Point (0.,0.,1.);
        nodes[26] = Point (0.,0.,0.);
        return;
      }
    case PRISM6:
      {
        nodes.resize(6);
        nodes[0] = Point (0.,0.,-1.);
        nodes[1] = Point (1.,0.,-1.);
        nodes[2] = Point (0.,1.,-1.);
        nodes[3] = Point (0.,0.,1.);
        nodes[4] = Point (1.,0.,1.);
        nodes[5] = Point (0.,1.,1.);
        return;
      }
    case PRISM15:
      {
        nodes.resize(15);
        nodes[0] = Point (0.,0.,-1.);
        nodes[1] = Point (1.,0.,-1.);
        nodes[2] = Point (0.,1.,-1.);
        nodes[3] = Point (0.,0.,1.);
        nodes[4] = Point (1.,0.,1.);
        nodes[5] = Point (0.,1.,1.);
        nodes[6] = Point (.5,0.,-1.);
        nodes[7] = Point (.5,.5,-1.);
        nodes[8] = Point (0.,.5,-1.);
        nodes[9] = Point (0.,0.,0.);
        nodes[10] = Point (1.,0.,0.);
        nodes[11] = Point (0.,1.,0.);
        nodes[12] = Point (.5,0.,1.);
        nodes[13] = Point (.5,.5,1.);
        nodes[14] = Point (0.,.5,1.);
        return;
      }
    case PRISM18:
      {
        nodes.resize(18);
        nodes[0] = Point (0.,0.,-1.);
        nodes[1] = Point (1.,0.,-1.);
        nodes[2] = Point (0.,1.,-1.);
        nodes[3] = Point (0.,0.,1.);
        nodes[4] = Point (1.,0.,1.);
        nodes[5] = Point (0.,1.,1.);
        nodes[6] = Point (.5,0.,-1.);
        nodes[7] = Point (.5,.5,-1.);
        nodes[8] = Point (0.,.5,-1.);
        nodes[9] = Point (0.,0.,0.);
        nodes[10] = Point (1.,0.,0.);
        nodes[11] = Point (0.,1.,0.);
        nodes[12] = Point (.5,0.,1.);
        nodes[13] = Point (.5,.5,1.);
        nodes[14] = Point (0.,.5,1.);
        nodes[15] = Point (.5,0.,0.);
        nodes[16] = Point (.5,.5,0.);
        nodes[17] = Point (0.,.5,0.);
        return;
      }
    case PYRAMID5:
      {
        nodes.resize(5);
        nodes[0] = Point (-1.,-1.,0.);
        nodes[1] = Point (1.,-1.,0.);
        nodes[2] = Point (1.,1.,0.);
        nodes[3] = Point (-1.,1.,0.);
        nodes[4] = Point (0.,0.,1.);
        return;
      }
    case PYRAMID13:
      {
        nodes.resize(13);

        // base corners
        nodes[0] = Point (-1.,-1.,0.);
        nodes[1] = Point (1.,-1.,0.);
        nodes[2] = Point (1.,1.,0.);
        nodes[3] = Point (-1.,1.,0.);

        // apex
        nodes[4] = Point (0.,0.,1.);

        // base midedge
        nodes[5] = Point (0.,-1.,0.);
        nodes[6] = Point (1.,0.,0.);
        nodes[7] = Point (0.,1.,0.);
        nodes[8] = Point (-1,0.,0.);

        // lateral midedge
        nodes[9] = Point (-.5,-.5,.5);
        nodes[10] = Point (.5,-.5,.5);
        nodes[11] = Point (.5,.5,.5);
        nodes[12] = Point (-.5,.5,.5);

        return;
      }
    case PYRAMID14:
      {
        nodes.resize(14);

        // base corners
        nodes[0] = Point (-1.,-1.,0.);
        nodes[1] = Point (1.,-1.,0.);
        nodes[2] = Point (1.,1.,0.);
        nodes[3] = Point (-1.,1.,0.);

        // apex
        nodes[4] = Point (0.,0.,1.);

        // base midedge
        nodes[5] = Point (0.,-1.,0.);
        nodes[6] = Point (1.,0.,0.);
        nodes[7] = Point (0.,1.,0.);
        nodes[8] = Point (-1,0.,0.);

        // lateral midedge
        nodes[9] = Point (-.5,-.5,.5);
        nodes[10] = Point (.5,-.5,.5);
        nodes[11] = Point (.5,.5,.5);
        nodes[12] = Point (-.5,.5,.5);

        // base center
        nodes[13] = Point (0.,0.,0.);

        return;
      }

    default:
      libmesh_error_msg("ERROR: Unknown element type " << itemType);
    }
}
const std::vector<std::vector<Point> >& libMesh::FEAbstract::get_tangents ( ) const [inline]
Returns:
the tangent vectors for face integration.

Definition at line 368 of file fe_abstract.h.

References _fe_map.

  { return this->_fe_map->get_tangents(); }
Returns:
the element type that the current shape functions have been calculated for. Useful in determining when shape functions must be recomputed.

Definition at line 407 of file fe_abstract.h.

References elem_type.

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

Definition at line 227 of file fe_abstract.h.

References _fe_map.

Referenced by libMesh::ExactSolution::_compute_error(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), and libMesh::ProjectFEMSolution::operator()().

  { return this->_fe_map->get_xyz(); }
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name) [inline, protected, inherited]

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 163 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.first++;
}
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name) [inline, protected, inherited]

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 176 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.second++;
}
virtual bool libMesh::FEAbstract::is_hierarchic ( ) const [pure virtual]
Returns:
true if the finite element's higher order shape functions are hierarchic

Implemented in libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, libMesh::FE< Dim, LAGRANGE_VEC >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, and libMesh::FE< Dim, T >.

static unsigned int libMesh::ReferenceCounter::n_objects ( ) [inline, static, inherited]

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 79 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

Referenced by libMesh::LibMeshInit::~LibMeshInit().

  { return _n_objects; }
bool libMesh::FEAbstract::on_reference_element ( const Point p,
const ElemType  t,
const Real  eps = TOLERANCE 
) [static]
Returns:
true if the point p is located on the reference element for element type t, false otherwise. Since we are doing floating point comparisons here the parameter eps can be specified to indicate a tolerance. For example, $ x \le 1 $ becomes $ x \le 1 + \epsilon $.

Definition at line 554 of file fe_abstract.C.

References libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::INFHEX8, libMesh::INFPRISM6, libMesh::NODEELEM, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::Real, libMesh::TET10, libMesh::TET4, libMesh::TRI3, and libMesh::TRI6.

Referenced by libMesh::FE< Dim, T >::inverse_map().

{
  libmesh_assert_greater_equal (eps, 0.);

  const Real xi   = p(0);
#if LIBMESH_DIM > 1
  const Real eta  = p(1);
#else
  const Real eta  = 0.;
#endif
#if LIBMESH_DIM > 2
  const Real zeta = p(2);
#else
  const Real zeta  = 0.;
#endif

  switch (t)
    {
    case NODEELEM:
      {
        return (!xi && !eta && !zeta);
      }
    case EDGE2:
    case EDGE3:
    case EDGE4:
      {
        // The reference 1D element is [-1,1].
        if ((xi >= -1.-eps) &&
            (xi <=  1.+eps))
          return true;

        return false;
      }


    case TRI3:
    case TRI6:
      {
        // The reference triangle is isocoles
        // and is bound by xi=0, eta=0, and xi+eta=1.
        if ((xi  >= 0.-eps) &&
            (eta >= 0.-eps) &&
            ((xi + eta) <= 1.+eps))
          return true;

        return false;
      }


    case QUAD4:
    case QUAD8:
    case QUAD9:
      {
        // The reference quadrilateral element is [-1,1]^2.
        if ((xi  >= -1.-eps) &&
            (xi  <=  1.+eps) &&
            (eta >= -1.-eps) &&
            (eta <=  1.+eps))
          return true;

        return false;
      }


    case TET4:
    case TET10:
      {
        // The reference tetrahedral is isocoles
        // and is bound by xi=0, eta=0, zeta=0,
        // and xi+eta+zeta=1.
        if ((xi   >= 0.-eps) &&
            (eta  >= 0.-eps) &&
            (zeta >= 0.-eps) &&
            ((xi + eta + zeta) <= 1.+eps))
          return true;

        return false;
      }


    case HEX8:
    case HEX20:
    case HEX27:
      {
        /*
          if ((xi   >= -1.) &&
          (xi   <=  1.) &&
          (eta  >= -1.) &&
          (eta  <=  1.) &&
          (zeta >= -1.) &&
          (zeta <=  1.))
          return true;
        */

        // The reference hexahedral element is [-1,1]^3.
        if ((xi   >= -1.-eps) &&
            (xi   <=  1.+eps) &&
            (eta  >= -1.-eps) &&
            (eta  <=  1.+eps) &&
            (zeta >= -1.-eps) &&
            (zeta <=  1.+eps))
          {
            //    libMesh::out << "Strange Point:\n";
            //    p.print();
            return true;
          }

        return false;
      }

    case PRISM6:
    case PRISM15:
    case PRISM18:
      {
        // Figure this one out...
        // inside the reference triange with zeta in [-1,1]
        if ((xi   >=  0.-eps) &&
            (eta  >=  0.-eps) &&
            (zeta >= -1.-eps) &&
            (zeta <=  1.+eps) &&
            ((xi + eta) <= 1.+eps))
          return true;

        return false;
      }


    case PYRAMID5:
    case PYRAMID13:
    case PYRAMID14:
      {
        // Check that the point is on the same side of all the faces
        // by testing whether:
        //
        // n_i.(x - x_i) <= 0
        //
        // for each i, where:
        //   n_i is the outward normal of face i,
        //   x_i is a point on face i.
        if ((-eta - 1. + zeta <= 0.+eps) &&
            (  xi - 1. + zeta <= 0.+eps) &&
            ( eta - 1. + zeta <= 0.+eps) &&
            ( -xi - 1. + zeta <= 0.+eps) &&
            (            zeta >= 0.-eps))
          return true;

        return false;
      }

#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
    case INFHEX8:
      {
        // The reference infhex8 is a [-1,1]^3.
        if ((xi   >= -1.-eps) &&
            (xi   <=  1.+eps) &&
            (eta  >= -1.-eps) &&
            (eta  <=  1.+eps) &&
            (zeta >= -1.-eps) &&
            (zeta <=  1.+eps))
          {
            return true;
          }
        return false;
      }

    case INFPRISM6:
      {
        // inside the reference triange with zeta in [-1,1]
        if ((xi   >=  0.-eps) &&
            (eta  >=  0.-eps) &&
            (zeta >= -1.-eps) &&
            (zeta <=  1.+eps) &&
            ((xi + eta) <= 1.+eps))
          {
            return true;
          }

        return false;
      }
#endif

    default:
      libmesh_error_msg("ERROR: Unknown element type " << t);
    }

  // If we get here then the point is _not_ in the
  // reference element.   Better return false.

  return false;
}
virtual void libMesh::FEAbstract::print_d2phi ( std::ostream &  os) const [pure virtual]

Prints the value of each shape function's second derivatives at each quadrature point. Implement in derived class since this depends on whether the element is vector-valued or not.

Implemented in libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.

virtual void libMesh::FEAbstract::print_dphi ( std::ostream &  os) const [pure virtual]

Prints the value of each shape function's derivative at each quadrature point. Implement in derived class since this depends on whether the element is vector-valued or not.

Implemented in libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.

Referenced by print_info().

void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out) [static, inherited]

Prints the reference information, by default to libMesh::out.

Definition at line 88 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

Referenced by libMesh::LibMeshInit::~LibMeshInit().

void libMesh::FEAbstract::print_info ( std::ostream &  os) const

Prints all the relevant information about the current element.

Definition at line 761 of file fe_abstract.C.

References print_dphi(), print_JxW(), print_phi(), and print_xyz().

Referenced by libMesh::operator<<().

{
  os << "phi[i][j]: Shape function i at quadrature pt. j" << std::endl;
  this->print_phi(os);

  os << "dphi[i][j]: Shape function i's gradient at quadrature pt. j" << std::endl;
  this->print_dphi(os);

  os << "XYZ locations of the quadrature pts." << std::endl;
  this->print_xyz(os);

  os << "Values of JxW at the quadrature pts." << std::endl;
  this->print_JxW(os);
}
void libMesh::FEAbstract::print_JxW ( std::ostream &  os) const

Prints the Jacobian times the weight for each quadrature point.

Definition at line 748 of file fe_abstract.C.

References _fe_map.

Referenced by print_info().

{
  this->_fe_map->print_JxW(os);
}
virtual void libMesh::FEAbstract::print_phi ( std::ostream &  os) const [pure virtual]

Prints the value of each shape function at each quadrature point. Implement in derived class since this depends on whether the element is vector-valued or not.

Implemented in libMesh::FEGenericBase< OutputType >, and libMesh::FEGenericBase< FEOutputType< T >::type >.

Referenced by print_info().

void libMesh::FEAbstract::print_xyz ( std::ostream &  os) const

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

Definition at line 755 of file fe_abstract.C.

References _fe_map.

Referenced by print_info().

{
  this->_fe_map->print_xyz(os);
}
virtual void libMesh::FEAbstract::reinit ( const Elem elem,
const std::vector< Point > *const  pts = NULL,
const std::vector< Real > *const  weights = NULL 
) [pure virtual]

This is at the core of this class. Use this for each new element in the mesh. Reinitializes the requested physical element-dependent data based on the current element elem. By default the element data are computed at the quadrature points specified by the quadrature rule qrule, but any set of points on the reference element may be specified in the optional argument pts.

Note that the FE classes decide which data to initialize based on which accessor functions such as get_phi() or get_d2phi() have been called, so all such accessors should be called before the first reinit().

Implemented in libMesh::FEXYZ< Dim >, libMesh::FESubdivision, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, and libMesh::FE< Dim, LAGRANGE_VEC >.

Referenced by libMesh::ExactSolution::_compute_error(), libMesh::FEMContext::build_new_fe(), and libMesh::ExactErrorEstimator::find_squared_element_error().

virtual void libMesh::FEAbstract::reinit ( const Elem elem,
const unsigned int  side,
const Real  tolerance = TOLERANCE,
const std::vector< Point > *const  pts = NULL,
const std::vector< Real > *const  weights = NULL 
) [pure virtual]

Reinitializes all the physical element-dependent data based on the side of the element elem. The tolerance paremeter is passed to the involved call to inverse_map(). By default the element data are computed at the quadrature points specified by the quadrature rule qrule, but any set of points on the reference side element may be specified in the optional argument pts.

Implemented in libMesh::FEXYZ< Dim >, libMesh::FESubdivision, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, libMesh::FE< Dim, LAGRANGE_VEC >, libMesh::FEXYZ< Dim >, and libMesh::FEXYZ< Dim >.

virtual bool libMesh::FEAbstract::shapes_need_reinit ( ) const [protected, pure virtual]
Returns:
true when the shape functions (for this FEFamily) depend on the particular element, and therefore needs to be re-initialized for each new element. false otherwise.

Implemented in libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::InfFE< Dim, T_radial, T_map >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< 2, SUBDIVISION >, libMesh::FE< Dim, HIERARCHIC >, libMesh::FE< Dim, SCALAR >, libMesh::FE< Dim, L2_LAGRANGE >, libMesh::FE< Dim, NEDELEC_ONE >, libMesh::FE< Dim, HERMITE >, libMesh::FE< Dim, CLOUGH >, libMesh::FE< Dim, MONOMIAL >, libMesh::FE< Dim, XYZ >, libMesh::FE< Dim, LAGRANGE >, libMesh::FE< Dim, L2_HIERARCHIC >, libMesh::FE< Dim, LAGRANGE_VEC >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, libMesh::FE< Dim, T >, and libMesh::FE< Dim, T >.

virtual void libMesh::FEAbstract::side_map ( const Elem elem,
const Elem side,
const unsigned int  s,
const std::vector< Point > &  reference_side_points,
std::vector< Point > &  reference_points 
) [pure virtual]

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const FEAbstract fe 
) [friend]

Same as above, but allows you to print to a stream.

Definition at line 777 of file fe_abstract.C.

{
  fe.print_info(os);
  return os;
}

Member Data Documentation

bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited]

Flag to control whether reference count information is printed when print_info is called.

Definition at line 137 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 131 of file reference_counter.h.

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited]

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 126 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

unsigned int libMesh::FEAbstract::_p_level [protected]

The p refinement level the current data structures are set up for.

Definition at line 570 of file fe_abstract.h.

Referenced by get_order(), and get_p_level().

bool libMesh::FEAbstract::calculate_curl_phi [mutable, protected]

Should we calculate shape function curls?

Definition at line 541 of file fe_abstract.h.

Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_curl_phi().

bool libMesh::FEAbstract::calculate_div_phi [mutable, protected]

Should we calculate shape function divergences?

Definition at line 546 of file fe_abstract.h.

Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_div_phi().

bool libMesh::FEAbstract::calculate_dphiref [mutable, protected]

Should we calculate reference shape function gradients?

Definition at line 551 of file fe_abstract.h.

Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_curl_phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phideta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidetadzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidx2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxi2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxidzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidy2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidydz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidz2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidzeta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_div_phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidx(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidxi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidz(), and libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidzeta().

bool libMesh::FEAbstract::calculate_phi [mutable, protected]

Should we calculate shape functions?

Definition at line 526 of file fe_abstract.h.

Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_phi(), and libMesh::FESubdivision::init_shape_functions().

bool libMesh::FEAbstract::calculations_started [mutable, protected]

Have calculations with this object already been started? Then all get_* functions should already have been called.

Definition at line 521 of file fe_abstract.h.

Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::get_curl_phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phideta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidetadzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidx2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxdz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxi2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidxidzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidy2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidydz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidz2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_d2phidzeta2(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_div_phi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphideta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidx(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidxi(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidy(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidz(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_dphidzeta(), libMesh::FEGenericBase< FEOutputType< T >::type >::get_phi(), and libMesh::FESubdivision::init_shape_functions().

const unsigned int libMesh::FEAbstract::dim [protected]

The dimensionality of the object

Definition at line 515 of file fe_abstract.h.

The element type the current data structures are set up for.

Definition at line 564 of file fe_abstract.h.

Referenced by libMesh::FESubdivision::attach_quadrature_rule(), and get_type().

The finite element type for this object. Note that this should be constant for the object.

Definition at line 558 of file fe_abstract.h.

Referenced by compute_node_constraints(), compute_periodic_node_constraints(), get_family(), get_fe_type(), get_order(), libMesh::InfFE< Dim, T_radial, T_map >::InfFE(), and libMesh::FESubdivision::init_shape_functions().

A pointer to the quadrature rule employed

Definition at line 575 of file fe_abstract.h.

Referenced by libMesh::FESubdivision::attach_quadrature_rule().

A flag indicating if current data structures correspond to quadrature rule points

Definition at line 581 of file fe_abstract.h.


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