| 
    Intrepid
    
   | 
 
Implementation of an H(grad)-compatible FEM basis of degree 2 on a Pyramid cell. More...
#include <Intrepid_HGRAD_PYR_I2_FEM.hpp>
  
 Public Member Functions | |
| Basis_HGRAD_PYR_I2_FEM () | |
| Constructor.  | |
| void | getValues (ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const | 
| FEM basis evaluation on a reference Pyramid cell.  More... | |
| void | getValues (ArrayScalar &outputValues, const ArrayScalar &inputPoints, const ArrayScalar &cellVertices, const EOperator operatorType=OPERATOR_VALUE) const | 
| FVD basis evaluation: invocation of this method throws an exception.  | |
  Public Member Functions inherited from Intrepid::Basis< Scalar, ArrayScalar > | |
| virtual | ~Basis () | 
| Destructor.  | |
| virtual int | getCardinality () const | 
| Returns cardinality of the basis.  More... | |
| virtual int | getDegree () const | 
| Returns the degree of the basis.  More... | |
| virtual const shards::CellTopology | getBaseCellTopology () const | 
| Returns the base cell topology for which the basis is defined. See Shards documentation http://trilinos.sandia.gov/packages/shards for definition of base cell topology.  More... | |
| virtual EBasis | getBasisType () const | 
| Returns the basis type.  More... | |
| virtual ECoordinates | getCoordinateSystem () const | 
| Returns the type of coordinate system for which the basis is defined.  More... | |
| virtual int | getDofOrdinal (const int subcDim, const int subcOrd, const int subcDofOrd) | 
| DoF tag to ordinal lookup.  More... | |
| 
virtual const std::vector < std::vector< std::vector < int > > > &  | getDofOrdinalData () | 
| DoF tag to ordinal data structure.  | |
| virtual const std::vector< int > & | getDofTag (const int dofOrd) | 
| DoF ordinal to DoF tag lookup.  More... | |
| virtual const std::vector < std::vector< int > > &  | getAllDofTags () | 
| Retrieves all DoF tags.  More... | |
Private Member Functions | |
| void | initializeTags () | 
| Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.  | |
Additional Inherited Members | |
  Protected Attributes inherited from Intrepid::Basis< Scalar, ArrayScalar > | |
| int | basisCardinality_ | 
| Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom.  | |
| int | basisDegree_ | 
| Degree of the largest complete polynomial space that can be represented by the basis.  | |
| shards::CellTopology | basisCellTopology_ | 
| Base topology of the cells for which the basis is defined. See the Shards package http://trilinos.sandia.gov/packages/shards for definition of base cell topology.  | |
| EBasis | basisType_ | 
| Type of the basis.  | |
| ECoordinates | basisCoordinates_ | 
| The coordinate system for which the basis is defined.  | |
| bool | basisTagsAreSet_ | 
| "true" if tagToOrdinal_ and ordinalToTag_ have been initialized  | |
| std::vector< std::vector< int > > | ordinalToTag_ | 
| DoF ordinal to tag lookup table.  More... | |
| std::vector< std::vector < std::vector< int > > >  | tagToOrdinal_ | 
| DoF tag to ordinal lookup table.  More... | |
Implementation of an H(grad)-compatible FEM basis of degree 2 on a Pyramid cell.
      Implements Lagrangian basis of degree 2 on the reference Pyramid cell. The basis has
      cardinality 13 and spans an INCOMPLETE bi-quadratic polynomial space. Basis functions are dual 
      to a unisolvent set of degrees-of-freedom (DoF) defined and enumerated as follows:
================================================================================================= | | degree-of-freedom-tag table | | | DoF |----------------------------------------------------------| DoF definition | | ordinal | subc dim | subc ordinal | subc DoF ord |subc num DoF | | |=========|==============|==============|==============|=============|===========================| | 0 | 0 | 0 | 0 | 1 | L_0(u) = u(-1,-1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 1 | 0 | 1 | 0 | 1 | L_1(u) = u( 1,-1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 2 | 0 | 2 | 0 | 1 | L_2(u) = u( 1, 1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 3 | 0 | 3 | 0 | 1 | L_3(u) = u(-1, 1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 4 | 0 | 4 | 0 | 1 | L_4(u) = u( 0, 0, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| |---------|--------------|--------------|--------------|-------------|---------------------------| | 5 | 1 | 0 | 0 | 1 | L_5(u) = u( 0,-1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 6 | 1 | 1 | 0 | 1 | L_6(u) = u( 1, 0, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 7 | 1 | 2 | 0 | 1 | L_7(u) = u( 0, 1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 8 | 1 | 3 | 0 | 1 | L_8(u) = u(-1, 0, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 9 | 1 | 4 | 0 | 1 | L_9(u) = u(-1/2,-1/2,1/2) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 10 | 1 | 5 | 0 | 1 | L_10(u)= u( 1/2,-1/2,1/2) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 11 | 1 | 6 | 0 | 1 | L_11(u)= u( 1/2, 1/2,1/2) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 12 | 1 | 7 | 0 | 1 | L_12(u)= u(-1/2, 1/2,1/2) | |=========|==============|==============|==============|=============|===========================| | MAX | maxScDim=1 | maxScOrd=7 | maxDfOrd=0 | - | | |=========|==============|==============|==============|=============|===========================|
\remark   Ordering of DoFs follows the node order in Pyramid<13> topology. 
\remark  WARNING: Does not satisfy a patch test for arbitrarily deformed pyramids, but only
                  for pyramids with bottom quadrilateral nodes forming a plane. 
Definition at line 105 of file Intrepid_HGRAD_PYR_I2_FEM.hpp.
      
  | 
  virtual | 
FEM basis evaluation on a reference Pyramid cell.
    Returns values of <var>operatorType</var> acting on FEM basis functions for a set of
    points in the <strong>reference Pyramid</strong> cell. For rank and dimensions of 
    I/O array arguments see Section \ref basis_md_array_sec .
| outputValues | [out] - rank-2 or 3 array with the computed basis values | 
| inputPoints | [in] - rank-2 array with dimensions (P,D) containing reference points | 
| operatorType | [in] - operator applied to basis functions | 
Implements Intrepid::Basis< Scalar, ArrayScalar >.
Definition at line 104 of file Intrepid_HGRAD_PYR_I2_FEMDef.hpp.
 1.8.5