Intrepid
Public Member Functions | Protected Attributes | Private Member Functions | List of all members
Intrepid::Basis< Scalar, ArrayScalar > Class Template Referenceabstract

An abstract base class that defines interface for concrete basis implementations for Finite Element (FEM) and Finite Volume/Finite Difference (FVD) discrete spaces. More...

#include <Intrepid_Basis.hpp>

Inheritance diagram for Intrepid::Basis< Scalar, ArrayScalar >:
Intrepid::Basis_HCURL_HEX_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HCURL_QUAD_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HCURL_TET_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HCURL_TET_In_FEM< Scalar, ArrayScalar > Intrepid::Basis_HCURL_TRI_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HCURL_TRI_In_FEM< Scalar, ArrayScalar > Intrepid::Basis_HCURL_WEDGE_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HDIV_HEX_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HDIV_QUAD_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HDIV_TET_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HDIV_TET_In_FEM< Scalar, ArrayScalar > Intrepid::Basis_HDIV_TRI_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HDIV_TRI_In_FEM< Scalar, ArrayScalar > Intrepid::Basis_HDIV_WEDGE_I1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_HEX_C1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_HEX_C2_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_HEX_I2_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_LINE_C1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_LINE_Cn_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_LINE_Hermite_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_POLY_C1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_PYR_C1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_PYR_I2_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_QUAD_C1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_QUAD_C2_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TET_C1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TET_C2_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TET_Cn_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TET_Cn_FEM_ORTH< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TET_COMP12_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TRI_C1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TRI_C2_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TRI_Cn_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_TRI_Cn_FEM_ORTH< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_WEDGE_C1_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_WEDGE_C2_FEM< Scalar, ArrayScalar > Intrepid::Basis_HGRAD_WEDGE_I2_FEM< Scalar, ArrayScalar > Intrepid::TensorBasis< Scalar, ArrayScalar >

Public Member Functions

virtual ~Basis ()
 Destructor.
 
virtual void getValues (ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const =0
 Evaluation of a FEM basis on a reference cell. More...
 
virtual void getValues (ArrayScalar &outputValues, const ArrayScalar &inputPoints, const ArrayScalar &cellVertices, const EOperator operatorType=OPERATOR_VALUE) const =0
 Evaluation of an FVD basis evaluation on a physical cell. More...
 
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...
 

Protected Attributes

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...
 

Private Member Functions

virtual void initializeTags ()=0
 Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.
 

Detailed Description

template<class Scalar, class ArrayScalar>
class Intrepid::Basis< Scalar, ArrayScalar >

An abstract base class that defines interface for concrete basis implementations for Finite Element (FEM) and Finite Volume/Finite Difference (FVD) discrete spaces.

      A FEM basis spans a discrete space whose type can be either COMPLETE or INCOMPLETE.  
      FEM basis functions are always defined on a reference cell and are dual to a unisolvent 
      set of degrees-of-freedom (DoF). FEM basis requires cell topology with a reference cell.

      An FVD basis spans a discrete space whose type is typically BROKEN. The basis functions 
      are defined directly on the physical cell and are dual to a set of DoFs on that cell. 
      As a result, FVD bases require the vertex coordinates of the physical cell but the cell
      itself is not required to have a reference cell.

      Every DoF and its corresponding basis function from a given FEM or FVD basis set is 
      assigned an ordinal number which which specifies its numerical position in the DoF set, 
      and a 4-field DoF tag whose first 3 fields establish association between the DoF and a 
      subcell of particular dimension, and the last field gives the total number of basis
      functions associated with that subcell; see Section \ref basis_dof_tag_ord_sec for details.
Remarks
To limit memory use by factory-type objects (basis factories will be included in future releases of Intrepid), tag data is not initialized by basis ctors, instead, whenever a function that requires tag data is invoked for a first time, it calls initializeTags() to fill ordinalToTag_ and tagToOrdinal_. Because tag data is basis specific, every concrete basis class requires its own implementation of initializeTags().
Todo:
restore test for inclusion of reference points in their resective reference cells in getValues_HGRAD_Args, getValues_CURL_Args, getValues_DIV_Args

Definition at line 89 of file Intrepid_Basis.hpp.

Member Function Documentation

template<class Scalar, class ArrayScalar>
const std::vector< std::vector< int > > & Basis::getAllDofTags ( )
virtual

Retrieves all DoF tags.

Returns
reference to a vector of vectors with dimensions (basisCardinality_, 4) such that
  • element [DofOrd][0] = tag field 0 for the DoF with the specified ordinal
  • element [DofOrd][1] = tag field 1 for the DoF with the specified ordinal
  • element [DofOrd][2] = tag field 2 for the DoF with the specified ordinal
  • element [DofOrd][3] = tag field 3 for the DoF with the specified ordinal

Definition at line 89 of file Intrepid_BasisDef.hpp.

Referenced by Intrepid::Basis_HGRAD_QUAD_Cn_FEM< Scalar, ArrayScalar >::initializeTags(), and Intrepid::Basis_HGRAD_HEX_Cn_FEM< Scalar, ArrayScalar >::initializeTags().

template<class Scalar, class ArrayScalar>
const shards::CellTopology Basis::getBaseCellTopology ( ) const
inlinevirtual
template<class Scalar, class ArrayScalar>
EBasis Basis::getBasisType ( ) const
inlinevirtual

Returns the basis type.

Returns
Basis type

Definition at line 106 of file Intrepid_BasisDef.hpp.

template<class Scalar, class ArrayScalar>
int Basis::getCardinality ( ) const
inlinevirtual
template<class Scalar, class ArrayScalar>
ECoordinates Basis::getCoordinateSystem ( ) const
inlinevirtual

Returns the type of coordinate system for which the basis is defined.

Returns
Type of the coordinate system (Cartesian, polar, R-Z, etc.).

Definition at line 124 of file Intrepid_BasisDef.hpp.

template<class Scalar, class ArrayScalar>
int Basis::getDegree ( ) const
inlinevirtual

Returns the degree of the basis.

Returns
max. degree of the complete polynomials that can be represented by the basis.

Definition at line 118 of file Intrepid_BasisDef.hpp.

template<class Scalar, class ArrayScalar>
int Basis::getDofOrdinal ( const int  subcDim,
const int  subcOrd,
const int  subcDofOrd 
)
virtual

DoF tag to ordinal lookup.

Parameters
subcDim[in] - tag field 0: dimension of the subcell associated with the DoF
subcOrd[in] - tag field 1: ordinal of the subcell defined by cell topology
subcDofOrd[in] - tag field 2: ordinal of the DoF relative to the subcell.
Returns
the DoF ordinal corresponding to the specified DoF tag data.

Definition at line 51 of file Intrepid_BasisDef.hpp.

template<class Scalar, class ArrayScalar>
const std::vector< int > & Basis::getDofTag ( const int  dofOrd)
virtual

DoF ordinal to DoF tag lookup.

Parameters
dofOrd[in] - ordinal of the DoF whose tag is being retrieved
Returns
reference to a vector with dimension (4) such that
  • element [0] = tag field 0 -> dim. of the subcell associated with the specified DoF
  • element [1] = tag field 1 -> ordinal of the subcell defined by cell topology
  • element [2] = tag field 2 -> ordinal of the specified DoF relative to the subcell
  • element [3] = tag field 3 -> total number of DoFs associated with the subcell

Definition at line 78 of file Intrepid_BasisDef.hpp.

template<class Scalar, class ArrayScalar>
virtual void Intrepid::Basis< Scalar, ArrayScalar >::getValues ( ArrayScalar &  outputValues,
const ArrayScalar &  inputPoints,
const EOperator  operatorType 
) const
pure virtual

Evaluation of a FEM basis on a reference cell.

    Returns values of <var>operatorType</var> acting on FEM basis functions for a set of
    points in the <strong>reference cell</strong> for which the basis is defined. 
Parameters
outputValues[out] - variable rank array with the basis values
inputPoints[in] - rank-2 array (P,D) with the evaluation points
operatorType[in] - the operator acting on the basis functions
Remarks
For rank and dimension specifications of the output array see Section MD array template arguments for basis methods. Dimensions of ArrayScalar arguments are checked at runtime if HAVE_INTREPID_DEBUG is defined.
A FEM basis spans a COMPLETE or INCOMPLETE polynomial space on the reference cell which is a smooth function space. Thus, all operator types that are meaningful for the approximated function space are admissible. When the order of the operator exceeds the degree of the basis, the output array is filled with the appropriate number of zeros.

Implemented in Intrepid::Basis_HGRAD_HEX_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_Hermite_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, Intrepid::FieldContainer< Scalar > >, Intrepid::Basis_HGRAD_HEX_I2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_WEDGE_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_HEX_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_WEDGE_I2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_WEDGE_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_HEX_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_PYR_I2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_TET_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_QUAD_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_WEDGE_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_TET_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_TET_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_TRI_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_COMP12_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_QUAD_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_TET_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_TRI_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TRI_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_TRI_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_QUAD_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_QUAD_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_HEX_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_HEX_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_QUAD_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_HEX_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_HEX_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_TRI_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TRI_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_PYR_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_WEDGE_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_QUAD_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TRI_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_QUAD_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_Cn_FEM_ORTH< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_Cn_FEM_ORTH< Scalar, Intrepid::FieldContainer< Scalar > >, Intrepid::Basis_HGRAD_TRI_Cn_FEM_ORTH< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TRI_Cn_FEM_ORTH< Scalar, Intrepid::FieldContainer< Scalar > >, and Intrepid::Basis_HGRAD_POLY_C1_FEM< Scalar, ArrayScalar >.

Referenced by Intrepid::Basis_HGRAD_QUAD_Cn_FEM< Scalar, ArrayScalar >::getValues(), and Intrepid::Basis_HGRAD_HEX_Cn_FEM< Scalar, ArrayScalar >::getValues().

template<class Scalar, class ArrayScalar>
virtual void Intrepid::Basis< Scalar, ArrayScalar >::getValues ( ArrayScalar &  outputValues,
const ArrayScalar &  inputPoints,
const ArrayScalar &  cellVertices,
const EOperator  operatorType = OPERATOR_VALUE 
) const
pure virtual

Evaluation of an FVD basis evaluation on a physical cell.

    Returns values of <var>operatorType</var> acting on FVD basis functions for a set of
    points in the <strong>physical cell</strong> for which the FVD basis is defined. 
Parameters
outputValues[out] - variable rank array with the basis values
inputPoints[in] - rank-2 array (P,D) with the evaluation points
cellVertices[in] - rank-2 array (V,D) with the vertices of the physical cell
operatorType[in] - the operator acting on the basis functions
Remarks
For rank and dimension specifications of the output array see Section MD array template arguments for basis methods. Dimensions of ArrayScalar arguments are checked at runtime if HAVE_INTREPID_DEBUG is defined.
A typical FVD basis spans a BROKEN discrete space which is only piecewise smooth. For example, it could be a piecewise constant space defined with respect to a partition of the cell into simplices. Because differential operators are not meaningful for such spaces, the default operator type in this method is set to OPERATOR_VALUE.

Implemented in Intrepid::Basis_HGRAD_HEX_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_Hermite_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, Intrepid::FieldContainer< Scalar > >, Intrepid::Basis_HGRAD_HEX_I2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_WEDGE_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_HEX_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_WEDGE_I2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_WEDGE_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_HEX_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_PYR_I2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_TET_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_QUAD_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_WEDGE_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_TET_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_TET_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_TRI_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_COMP12_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_QUAD_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_TET_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_TRI_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TRI_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_TRI_I1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_QUAD_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_QUAD_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HDIV_HEX_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_HEX_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_QUAD_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_HEX_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_HEX_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HCURL_TRI_In_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TRI_C2_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_PYR_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_WEDGE_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_QUAD_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TRI_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_LINE_C1_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_QUAD_Cn_FEM< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_Cn_FEM_ORTH< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TET_Cn_FEM_ORTH< Scalar, Intrepid::FieldContainer< Scalar > >, Intrepid::Basis_HGRAD_TRI_Cn_FEM_ORTH< Scalar, ArrayScalar >, Intrepid::Basis_HGRAD_TRI_Cn_FEM_ORTH< Scalar, Intrepid::FieldContainer< Scalar > >, and Intrepid::Basis_HGRAD_POLY_C1_FEM< Scalar, ArrayScalar >.

Member Data Documentation

template<class Scalar, class ArrayScalar>
std::vector<std::vector<int> > Intrepid::Basis< Scalar, ArrayScalar >::ordinalToTag_
protected

DoF ordinal to tag lookup table.

    Rank-2 array with dimensions (basisCardinality_, 4) containing the DoF tags. This array 
    is left empty at instantiation and filled by initializeTags() only when tag data is
    requested. 
  • ordinalToTag_[DofOrd][0] = dim. of the subcell associated with the specified DoF
  • ordinalToTag_[DofOrd][1] = ordinal of the subcell defined in the cell topology
  • ordinalToTag_[DodOrd][2] = ordinal of the specified DoF relative to the subcell
  • ordinalToTag_[DofOrd][3] = total number of DoFs associated with the subcell

Definition at line 134 of file Intrepid_Basis.hpp.

template<class Scalar, class ArrayScalar>
std::vector<std::vector<std::vector<int> > > Intrepid::Basis< Scalar, ArrayScalar >::tagToOrdinal_
protected

DoF tag to ordinal lookup table.

    Rank-3 array with dimensions (maxScDim + 1, maxScOrd + 1, maxDfOrd + 1), i.e., the 
    columnwise maximums of the 1st three columns in the DoF tag table for the basis plus 1. 
    For every triple (subscDim, subcOrd, subcDofOrd) that is valid DoF tag data this array 
    stores the corresponding DoF ordinal. If the triple does not correspond to tag data, 
    the array stores -1. This array is left empty at instantiation and filled by 
    initializeTags() only when tag data is requested. 
  • tagToOrdinal_[subcDim][subcOrd][subcDofOrd] = Degree-of-freedom ordinal

Definition at line 147 of file Intrepid_Basis.hpp.


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