Intrepid2
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Intrepid2::Basis< DeviceType, OutputType, PointType > Class Template Reference

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

Public Types

using DeviceType = Device
 (Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_device (may be Kokkos::Serial, for example).
 
using ExecutionSpace = typename DeviceType::execution_space
 (Kokkos) Execution space for basis.
 
using OutputValueType = outputValueType
 Output value type for basis; default is double.
 
using PointValueType = pointValueType
 Point value type for basis; default is double.
 
using OrdinalViewType = Kokkos::View< ordinal_type, DeviceType >
 View type for ordinal.
 
using EBasisViewType = Kokkos::View< EBasis, DeviceType >
 View for basis type.
 
using ECoordinatesViewType = Kokkos::View< ECoordinates, DeviceType >
 View for coordinate system type.
 
using OrdinalTypeArray1DHost = Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace >
 View type for 1d host array.
 
using OrdinalTypeArray2DHost = Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace >
 View type for 2d host array.
 
using OrdinalTypeArray3DHost = Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace >
 View type for 3d host array.
 
using OrdinalTypeArrayStride1DHost = Kokkos::View< ordinal_type *, Kokkos::LayoutStride, Kokkos::HostSpace >
 View type for 1d host array.
 
using OrdinalTypeArray1D = Kokkos::View< ordinal_type *, DeviceType >
 View type for 1d device array.
 
using OrdinalTypeArray2D = Kokkos::View< ordinal_type **, DeviceType >
 View type for 2d device array.
 
using OrdinalTypeArray3D = Kokkos::View< ordinal_type ***, DeviceType >
 View type for 3d device array.
 
using OrdinalTypeArrayStride1D = Kokkos::View< ordinal_type *, Kokkos::LayoutStride, DeviceType >
 View type for 1d device array.
 
typedef ScalarTraits
< pointValueType >
::scalar_type 
scalarType
 Scalar type for point values.
 
using OutputViewType = Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType >
 View type for basis value output.
 
using PointViewType = Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType >
 View type for input points.
 
using ScalarViewType = Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType >
 View type for scalars.
 

Public Member Functions

OutputValueType getDummyOutputValue ()
 Dummy array to receive input arguments.
 
PointValueType getDummyPointValue ()
 Dummy array to receive input arguments.
 
Kokkos::DynRankView
< OutputValueType, DeviceType
allocateOutputView (const int numPoints, const EOperator operatorType=OPERATOR_VALUE) const
 Allocate a View container suitable for passing to the getValues() variant that accepts Kokkos DynRankViews as arguments (as opposed to the Intrepid2 BasisValues and PointValues containers). More...
 
virtual BasisValues
< OutputValueType, DeviceType
allocateBasisValues (TensorPoints< PointValueType, DeviceType > points, const EOperator operatorType=OPERATOR_VALUE) const
 Allocate BasisValues container suitable for passing to the getValues() variant that takes a TensorPoints container as argument. More...
 
virtual void getValues (const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
 Evaluation of a FEM basis on a reference cell. More...
 
virtual void getValues (OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
 
virtual void getValues (BasisValues< OutputValueType, DeviceType > outputValues, const TensorPoints< PointValueType, DeviceType > inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
 Evaluation of a FEM basis on a reference cell, using point and output value containers that allow preservation of tensor-product structure. More...
 
virtual void getValues (OutputViewType, const PointViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
 Evaluation of an FVD basis evaluation on a physical cell. More...
 
virtual void getDofCoords (ScalarViewType) const
 Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
 
virtual void getDofCoeffs (ScalarViewType) const
 Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space spanned by the basis, := P(dofCoords(i)) dofCoeffs(i) are the nodal coefficients associated to basis function i. More...
 
OrdinalTypeArray1DHost getFieldOrdinalsForDegree (OrdinalTypeArray1DHost &degrees) const
 For hierarchical bases, returns the field ordinals that have at most the specified degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees. More...
 
OrdinalTypeArray1DHost getFieldOrdinalsForH1Degree (OrdinalTypeArray1DHost &degrees) const
 For hierarchical bases, returns the field ordinals that have at most the specified H^1 degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees. More...
 
std::vector< int > getFieldOrdinalsForDegree (std::vector< int > &degrees) const
 For hierarchical bases, returns the field ordinals that have at most the specified degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees. More...
 
std::vector< int > getFieldOrdinalsForH1Degree (std::vector< int > &degrees) const
 For hierarchical bases, returns the field ordinals that have at most the specified H^1 degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees. More...
 
OrdinalTypeArray1DHost getPolynomialDegreeOfField (int fieldOrdinal) const
 For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array. More...
 
OrdinalTypeArray1DHost getH1PolynomialDegreeOfField (int fieldOrdinal) const
 For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array. More...
 
std::vector< int > getPolynomialDegreeOfFieldAsVector (int fieldOrdinal) const
 For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array. More...
 
std::vector< int > getH1PolynomialDegreeOfFieldAsVector (int fieldOrdinal) const
 For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array. More...
 
int getPolynomialDegreeLength () const
 For hierarchical bases, returns the number of entries required to specify the polynomial degree of a basis function.
 
virtual const char * getName () const
 Returns basis name. More...
 
virtual bool requireOrientation () const
 True if orientation is required.
 
ordinal_type getCardinality () const
 Returns cardinality of the basis. More...
 
ordinal_type getDegree () const
 Returns the degree of the basis. More...
 
EFunctionSpace getFunctionSpace () const
 Returns the function space for the basis. More...
 
shards::CellTopology getBaseCellTopology () const
 Returns the base cell topology for which the basis is defined. See Shards documentation https://trilinos.org/packages/shards for definition of base cell topology. More...
 
EBasis getBasisType () const
 Returns the basis type. More...
 
ECoordinates getCoordinateSystem () const
 Returns the type of coordinate system for which the basis is defined. More...
 
ordinal_type getDofCount (const ordinal_type subcDim, const ordinal_type subcOrd) const
 DoF count for specified subcell. More...
 
ordinal_type getDofOrdinal (const ordinal_type subcDim, const ordinal_type subcOrd, const ordinal_type subcDofOrd) const
 DoF tag to ordinal lookup. More...
 
virtual int getNumTensorialExtrusions () const
 returns the number of tensorial extrusions relative to the cell topology returned by getBaseCellTopology(). Base class returns 0; overridden by TensorBasis.
 
const OrdinalTypeArray3DHost getAllDofOrdinal () const
 DoF tag to ordinal data structure.
 
const OrdinalTypeArrayStride1DHost getDofTag (const ordinal_type dofOrd) const
 DoF ordinal to DoF tag lookup. More...
 
const OrdinalTypeArray2DHost getAllDofTags () const
 Retrieves all DoF tags. More...
 
virtual BasisPtr< DeviceType,
OutputValueType,
PointValueType
getSubCellRefBasis (const ordinal_type subCellDim, const ordinal_type subCellOrd) const
 returns the basis associated to a subCell. More...
 
ordinal_type getDomainDimension () const
 Returns the spatial dimension of the domain of the basis; this is equal to getBaseCellTopology().getDimension() + getNumTensorialExtrusions(). More...
 
virtual HostBasisPtr
< OutputValueType,
PointValueType
getHostBasis () const
 Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_type, but is otherwise identical to this. More...
 

Protected Member Functions

template<typename OrdinalTypeView3D , typename OrdinalTypeView2D , typename OrdinalTypeView1D >
void setOrdinalTagData (OrdinalTypeView3D &tagToOrdinal, OrdinalTypeView2D &ordinalToTag, const OrdinalTypeView1D tags, const ordinal_type basisCard, const ordinal_type tagSize, const ordinal_type posScDim, const ordinal_type posScOrd, const ordinal_type posDfOrd)
 Fills ordinalToTag_ and tagToOrdinal_ by basis-specific tag data. More...
 

Protected Attributes

ordinal_type basisCardinality_
 Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom.
 
ordinal_type 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 for definition of base cell topology. For TensorBasis subclasses, by default this the cell topology that is extruded (i.e., it is a lower-dimensional CellTopology than the space on which the tensor basis is defined). This allows tensor bases to be defined in higher dimensions than shards::CellTopology supports. TensorBasis subclasses can opt to use an equivalent shards CellTopology for basisCellTopology_, as well as using Intrepid2's tagging for tensor bases in dimensions up to 3, by calling TensorBasis::setShardsTopologyAndTags().
 
EBasis basisType_
 Type of the basis.
 
ECoordinates basisCoordinates_
 The coordinate system for which the basis is defined.
 
EFunctionSpace functionSpace_ = FUNCTION_SPACE_MAX
 The function space in which the basis is defined.
 
OrdinalTypeArray2DHost ordinalToTag_
 "true" if tagToOrdinal_ and ordinalToTag_ have been initialized More...
 
OrdinalTypeArray3DHost tagToOrdinal_
 DoF tag to ordinal lookup table. More...
 
Kokkos::DynRankView
< scalarType, DeviceType
dofCoords_
 Coordinates of degrees-of-freedom for basis functions defined in physical space.
 
Kokkos::DynRankView
< scalarType, DeviceType
dofCoeffs_
 Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space spanned by the basis, := P(dofCoords_(i)) dofCoeffs_(i) are the nodal coefficients associated to basis functions i. More...
 
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
 Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now. The number of entries per degree of freedom in this table depends on the basis type. For hypercubes, this will be the spatial dimension. We have not yet determined what this will be for simplices beyond 1D; there are not yet hierarchical simplicial bases beyond 1D in Intrepid2. More...
 
OrdinalTypeArray2DHost fieldOrdinalH1PolynomialDegree_
 H^1 polynomial degree for each degree of freedom. Only defined for hierarchical bases right now. The number of entries per degree of freedom in this table depends on the basis type. For hypercubes, this will be the spatial dimension. We have not yet determined what this will be for simplices beyond 1D; there are not yet hierarchical simplicial bases beyond 1D in Intrepid2. More...
 

Detailed Description

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
class Intrepid2::Basis< DeviceType, OutputType, PointType >

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.

    Note the use of Kokkos::LayoutStride as the layout for various View definitions (including the
    argument for \ref getValues() ).  A method expecting a LayoutStride view can accept, thanks to
    some conversion mechanisms defined elsewhere, a LayoutLeft view (the default for Cuda views in
    Kokkos) or a LayoutRight view (the default on most other platforms).  This does introduce some
    additional complexity when Views need to be allocated for temporary storage; see the method
    \ref getMatchingViewWithLabel() provided in \ref Intrepid_Utils.hpp.
Remarks
To limit memory use by factory-type objects (basis factories will be included in future releases of Intrepid2), 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 71 of file Intrepid2_Basis.hpp.

Member Function Documentation

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual BasisValues<OutputValueType,DeviceType> Intrepid2::Basis< DeviceType, OutputType, PointType >::allocateBasisValues ( TensorPoints< PointValueType, DeviceType points,
const EOperator  operatorType = OPERATOR_VALUE 
) const
inlinevirtual

Allocate BasisValues container suitable for passing to the getValues() variant that takes a TensorPoints container as argument.

The default implementation employs a trivial tensor-product structure, for compatibility across all bases. Subclasses that have non-trivial tensor-product structure should override. The basic exact-sequence operators are supported (VALUE, GRAD, DIV, CURL), as are the Dn operators (OPERATOR_D1 through OPERATOR_D10).

Definition at line 388 of file Intrepid2_Basis.hpp.

template<typename Device , typename outputValueType , typename pointValueType >
Kokkos::DynRankView< outputValueType, Device > Intrepid2::Basis< Device, outputValueType, pointValueType >::allocateOutputView ( const int  numPoints,
const EOperator  operatorType = OPERATOR_VALUE 
) const

Allocate a View container suitable for passing to the getValues() variant that accepts Kokkos DynRankViews as arguments (as opposed to the Intrepid2 BasisValues and PointValues containers).

Note that only the basic exact-sequence operators are supported at the moment: VALUE, GRAD, DIV, CURL.

Definition at line 860 of file Intrepid2_BasisDef.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::allocateBasisValues().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
const OrdinalTypeArray2DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::getAllDofTags ( ) const
inline

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 937 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
shards::CellTopology Intrepid2::Basis< DeviceType, OutputType, PointType >::getBaseCellTopology ( ) const
inline

Returns the base cell topology for which the basis is defined. See Shards documentation https://trilinos.org/packages/shards for definition of base cell topology.

Returns
Base cell topology

Definition at line 812 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getDomainDimension().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
EBasis Intrepid2::Basis< DeviceType, OutputType, PointType >::getBasisType ( ) const
inline

Returns the basis type.

Returns
Basis type

Definition at line 822 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
ordinal_type Intrepid2::Basis< DeviceType, OutputType, PointType >::getCardinality ( ) const
inline

Returns cardinality of the basis.

Returns
the number of basis functions in the basis

Definition at line 783 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::ProjectionTools< DeviceType >::getHCurlBasisCoeffs(), and Intrepid2::ProjectionTools< DeviceType >::getHDivBasisCoeffs().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
ECoordinates Intrepid2::Basis< DeviceType, OutputType, PointType >::getCoordinateSystem ( ) const
inline

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 832 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
ordinal_type Intrepid2::Basis< DeviceType, OutputType, PointType >::getDegree ( ) const
inline

Returns the degree of the basis.

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

Definition at line 793 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual void Intrepid2::Basis< DeviceType, OutputType, PointType >::getDofCoeffs ( ScalarViewType  ) const
inlinevirtual

Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space spanned by the basis, := P(dofCoords(i)) dofCoeffs(i) are the nodal coefficients associated to basis function i.

Rank-1 array for scalar basis with dimension (cardinality) Rank-2 array for vector basis with dimensions (cardinality, cell dimension)

Reimplemented in Intrepid2::Basis_HGRAD_HEX_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_C1_FEM< DeviceType, outputValueType, pointValueType >, and Intrepid2::Basis_HGRAD_LINE_C2_FEM< DeviceType, outputValueType, pointValueType >.

Definition at line 539 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
ordinal_type Intrepid2::Basis< DeviceType, OutputType, PointType >::getDofCount ( const ordinal_type  subcDim,
const ordinal_type  subcOrd 
) const
inline

DoF count for specified subcell.

Parameters
subcDim[in] - tag field 0: dimension of the subcell associated with the DoFs
subcOrd[in] - tag field 1: ordinal of the subcell defined by cell topology
Returns
the number of DoFs associated with the specified subcell.

Definition at line 844 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::ProjectionStruct< DeviceType, ValueType >::createHDivProjectionStruct(), Intrepid2::ProjectionTools< DeviceType >::getHCurlBasisCoeffs(), and Intrepid2::ProjectionTools< DeviceType >::getHDivBasisCoeffs().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
ordinal_type Intrepid2::Basis< DeviceType, OutputType, PointType >::getDofOrdinal ( const ordinal_type  subcDim,
const ordinal_type  subcOrd,
const ordinal_type  subcDofOrd 
) const
inline

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 870 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
const OrdinalTypeArrayStride1DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::getDofTag ( const ordinal_type  dofOrd) const
inline

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 919 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getDofCount().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
ordinal_type Intrepid2::Basis< DeviceType, OutputType, PointType >::getDomainDimension ( ) const
inline

Returns the spatial dimension of the domain of the basis; this is equal to getBaseCellTopology().getDimension() + getNumTensorialExtrusions().

Returns
The spatial dimension of the domain.

Definition at line 965 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
OrdinalTypeArray1DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::getFieldOrdinalsForDegree ( OrdinalTypeArray1DHost degrees) const
inline

For hierarchical bases, returns the field ordinals that have at most the specified degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees.

Parameters
degrees[in] - 1D host ordinal array of length specified by getPolynomialDegreeLength(), indicating what the maximum degree in each dimension should be
Returns
a 1D host ordinal array containing the ordinals of matching basis functions

Definition at line 551 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getFieldOrdinalsForDegree().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
std::vector<int> Intrepid2::Basis< DeviceType, OutputType, PointType >::getFieldOrdinalsForDegree ( std::vector< int > &  degrees) const
inline

For hierarchical bases, returns the field ordinals that have at most the specified degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees.

This variant takes a std::vector of polynomial degrees and returns a std::vector of field ordinals. It calls the other variant, which uses Kokkos Views on the host.

Parameters
degrees[in] - std::vector<int> of length specified by getPolynomialDegreeLength(), indicating what the maximum degree in each dimension should be
Returns
a std::vector<int> containing the ordinals of matching basis functions

Definition at line 619 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
OrdinalTypeArray1DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::getFieldOrdinalsForH1Degree ( OrdinalTypeArray1DHost degrees) const
inline

For hierarchical bases, returns the field ordinals that have at most the specified H^1 degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees.

Parameters
degrees[in] - 1D host ordinal array of length specified by getPolynomialDegreeLength(), indicating what the maximum degree in each dimension should be
Returns
a 1D host ordinal array containing the ordinals of matching basis functions

Definition at line 583 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getFieldOrdinalsForH1Degree().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
std::vector<int> Intrepid2::Basis< DeviceType, OutputType, PointType >::getFieldOrdinalsForH1Degree ( std::vector< int > &  degrees) const
inline

For hierarchical bases, returns the field ordinals that have at most the specified H^1 degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees.

This variant takes a std::vector of polynomial degrees and returns a std::vector of field ordinals. It calls the other variant, which uses Kokkos Views on the host.

Parameters
degrees[in] - std::vector<int> of length specified by getPolynomialDegreeLength(), indicating what the maximum degree in each dimension should be
Returns
a std::vector<int> containing the ordinals of matching basis functions

Definition at line 647 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
EFunctionSpace Intrepid2::Basis< DeviceType, OutputType, PointType >::getFunctionSpace ( ) const
inline

Returns the function space for the basis.

Returns
the function space.

Definition at line 802 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
OrdinalTypeArray1DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::getH1PolynomialDegreeOfField ( int  fieldOrdinal) const
inline

For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array.

Parameters
fieldOrdinal[in] - ordinal of the basis function whose polynomial degree is requested.
Returns
a 1D host array of length matching getPolynomialDegreeLength(), with the H^1 polynomial degree of the basis function in each dimension.

Definition at line 693 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getH1PolynomialDegreeOfFieldAsVector().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
std::vector<int> Intrepid2::Basis< DeviceType, OutputType, PointType >::getH1PolynomialDegreeOfFieldAsVector ( int  fieldOrdinal) const
inline

For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array.

Parameters
fieldOrdinal[in] - ordinal of the basis function whose polynomial degree is requested.
Returns
a std::vector<int> of length matching getPolynomialDegreeLength(), with the polynomial degree of the basis function in each dimension.

Definition at line 737 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual HostBasisPtr<OutputValueType, PointValueType> Intrepid2::Basis< DeviceType, OutputType, PointType >::getHostBasis ( ) const
inlinevirtual

Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_type, but is otherwise identical to this.

Returns
Pointer to the new Basis object.

Reimplemented in Intrepid2::HierarchicalBasis_HDIV_PYR< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::IntegratedLegendreBasis_HGRAD_PYR< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::HierarchicalBasis_HCURL_TET< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::HierarchicalBasis_HDIV_TET< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::IntegratedLegendreBasis_HGRAD_TET< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::IntegratedLegendreBasis_HGRAD_TRI< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::HierarchicalBasis_HCURL_TRI< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::LegendreBasis_HVOL_PYR< DeviceType, OutputScalar, PointScalar >, Intrepid2::IntegratedLegendreBasis_HGRAD_LINE< DeviceType, OutputScalar, PointScalar, defineVertexFunctions, useMinusOneToOneReferenceElement >, Intrepid2::LegendreBasis_HVOL_TET< DeviceType, OutputScalar, PointScalar >, Intrepid2::LegendreBasis_HVOL_TRI< DeviceType, OutputScalar, PointScalar >, Intrepid2::Basis_HGRAD_HEX_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::LegendreBasis_HVOL_LINE< DeviceType, OutputScalar, PointScalar >, Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_QUAD_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_QUAD_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_PYR_I2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_WEDGE_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_PYR_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_C0_FEM< DeviceType, outputValueType, pointValueType >, and Intrepid2::HierarchicalBasis_HDIV_TRI< DeviceType, OutputScalar, PointScalar, useCGBasis >.

Definition at line 975 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual const char* Intrepid2::Basis< DeviceType, OutputType, PointType >::getName ( ) const
inlinevirtual

Returns basis name.

Returns
the name of the basis

Reimplemented in Intrepid2::HierarchicalBasis_HDIV_PYR< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::IntegratedLegendreBasis_HGRAD_PYR< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::HierarchicalBasis_HCURL_TET< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::HierarchicalBasis_HDIV_TET< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::IntegratedLegendreBasis_HGRAD_TET< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::IntegratedLegendreBasis_HGRAD_TRI< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::HierarchicalBasis_HCURL_TRI< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::LegendreBasis_HVOL_PYR< DeviceType, OutputScalar, PointScalar >, Intrepid2::IntegratedLegendreBasis_HGRAD_LINE< DeviceType, OutputScalar, PointScalar, defineVertexFunctions, useMinusOneToOneReferenceElement >, Intrepid2::LegendreBasis_HVOL_TET< DeviceType, OutputScalar, PointScalar >, Intrepid2::Basis_HGRAD_HEX_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::LegendreBasis_HVOL_TRI< DeviceType, OutputScalar, PointScalar >, Intrepid2::Basis_HDIV_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::LegendreBasis_HVOL_LINE< DeviceType, OutputScalar, PointScalar >, Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_COMP12_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_QUAD_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_QUAD_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_PYR_I2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_WEDGE_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_PYR_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_C0_FEM< DeviceType, outputValueType, pointValueType >, and Intrepid2::HierarchicalBasis_HDIV_TRI< DeviceType, OutputScalar, PointScalar, useCGBasis >.

Definition at line 766 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
OrdinalTypeArray1DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::getPolynomialDegreeOfField ( int  fieldOrdinal) const
inline

For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array.

Parameters
fieldOrdinal[in] - ordinal of the basis function whose polynomial degree is requested.
Returns
a 1D host array of length matching getPolynomialDegreeLength(), with the polynomial degree of the basis function in each dimension.

Definition at line 671 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getPolynomialDegreeOfFieldAsVector().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
std::vector<int> Intrepid2::Basis< DeviceType, OutputType, PointType >::getPolynomialDegreeOfFieldAsVector ( int  fieldOrdinal) const
inline

For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array.

Parameters
fieldOrdinal[in] - ordinal of the basis function whose polynomial degree is requested.
Returns
a std::vector<int> of length matching getPolynomialDegreeLength(), with the polynomial degree of the basis function in each dimension.

Definition at line 716 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual BasisPtr<DeviceType, OutputValueType, PointValueType> Intrepid2::Basis< DeviceType, OutputType, PointType >::getSubCellRefBasis ( const ordinal_type  subCellDim,
const ordinal_type  subCellOrd 
) const
inlinevirtual

returns the basis associated to a subCell.

HGRAD case: The bases of the subCell are the restriction to the subCell of the bases of the parent cell. HCURL case: The bases of the subCell are the restriction to the subCell of the bases of the parent cell, projected onto the manifold tangent to the subCell HDIV case: The bases of the subCell are the restriction to the subCell of the bases of the parent cell, projected along the normal of the subCell

This method is not supported by all bases (e.g. bases defined on a line and HVOL bases).

Parameters
[in]subCellDim- dimension of subCell
[in]subCellOrd- position of the subCell among of the subCells having the same dimension
Returns
pointer to the subCell basis of dimension subCellDim and position subCellOrd

Reimplemented in Intrepid2::HierarchicalBasis_HDIV_PYR< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::IntegratedLegendreBasis_HGRAD_PYR< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::HierarchicalBasis_HCURL_TET< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::HierarchicalBasis_HDIV_TET< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::IntegratedLegendreBasis_HGRAD_TET< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::IntegratedLegendreBasis_HGRAD_TRI< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::HierarchicalBasis_HCURL_TRI< DeviceType, OutputScalar, PointScalar, useCGBasis >, Intrepid2::LegendreBasis_HVOL_PYR< DeviceType, OutputScalar, PointScalar >, Intrepid2::LegendreBasis_HVOL_TRI< DeviceType, OutputScalar, PointScalar >, Intrepid2::Basis_HGRAD_HEX_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_QUAD_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_QUAD_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_PYR_I2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_WEDGE_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_PYR_C1_FEM< DeviceType, outputValueType, pointValueType >, and Intrepid2::Basis_HGRAD_TRI_C1_FEM< DeviceType, outputValueType, pointValueType >.

Definition at line 956 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual void Intrepid2::Basis< DeviceType, OutputType, PointType >::getValues ( const ExecutionSpace ,
OutputViewType  ,
const PointViewType  ,
const EOperator  = OPERATOR_VALUE 
) const
inlinevirtual

Evaluation of a FEM basis on a reference cell.

Returns values of operatorType acting on FEM basis functions for a set of points in the reference cell for which the basis is defined.

Parameters
space[in] - execution space instance
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_INTREPID2_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.

Reimplemented in Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_COMP12_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM< serendipity, DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_C1_FEM< DeviceType, outputValueType, pointValueType >, and Intrepid2::Basis_HGRAD_LINE_C2_FEM< DeviceType, outputValueType, pointValueType >.

Definition at line 436 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::ProjectionTools< DeviceType >::getHCurlBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getHDivBasisCoeffs(), and Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getValues().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual void Intrepid2::Basis< DeviceType, OutputType, PointType >::getValues ( OutputViewType  outputValues,
const PointViewType  inputPoints,
const EOperator  operatorType = OPERATOR_VALUE 
) const
inlinevirtual
template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual void Intrepid2::Basis< DeviceType, OutputType, PointType >::getValues ( BasisValues< OutputValueType, DeviceType outputValues,
const TensorPoints< PointValueType, DeviceType inputPoints,
const EOperator  operatorType = OPERATOR_VALUE 
) const
inlinevirtual

Evaluation of a FEM basis on a reference cell, using point and output value containers that allow preservation of tensor-product structure.

Returns values of operatorType acting on FEM basis functions for a set of points in the reference cell for which the basis is defined.

Parameters
outputValues[out] - variable rank array with the basis values. Should be allocated using Basis::allocateBasisValues().
inputPoints[in] - rank-2 array (P,D) with the evaluation points. This should be allocated using Cubature::allocateCubaturePoints() and filled using Cubature::getCubature().
operatorType[in] - the operator acting on the basis function

The default implementation does not take advantage of any tensor-product structure; subclasses with tensor-product support must override allocateBasisValues() and this getValues() method.

Definition at line 466 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
virtual void Intrepid2::Basis< DeviceType, OutputType, PointType >::getValues ( OutputViewType  ,
const PointViewType  ,
const PointViewType  ,
const EOperator  = OPERATOR_VALUE 
) const
inlinevirtual

Evaluation of an FVD basis evaluation on a physical cell.

Returns values of operatorType acting on FVD basis functions for a set of points in the physical cell 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_INTREPID2_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.

Definition at line 510 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
template<typename OrdinalTypeView3D , typename OrdinalTypeView2D , typename OrdinalTypeView1D >
void Intrepid2::Basis< DeviceType, OutputType, PointType >::setOrdinalTagData ( OrdinalTypeView3D &  tagToOrdinal,
OrdinalTypeView2D &  ordinalToTag,
const OrdinalTypeView1D  tags,
const ordinal_type  basisCard,
const ordinal_type  tagSize,
const ordinal_type  posScDim,
const ordinal_type  posScOrd,
const ordinal_type  posDfOrd 
)
inlineprotected

Fills ordinalToTag_ and tagToOrdinal_ by basis-specific tag data.

Parameters
tagToOrdinal[out] - Lookup table for the DoF's ordinal by its tag
ordinalToTag[out] - Lookup table for the DoF's tag by its ordinal
tags[in] - a set of basis-dependent tags in flat (rank-1) array format.
basisCard[in] - cardinality of the basis
tagSize[in] - number of fields in a DoF tag
posScDim[in] - position in the tag, counting from 0, of the subcell dim
posScOrd[in] - position in the tag, counting from 0, of the subcell ordinal
posDfOrd[in] - position in the tag, counting from 0, of DoF ordinal relative to the subcell

Definition at line 266 of file Intrepid2_Basis.hpp.

Member Data Documentation

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
Kokkos::DynRankView<scalarType,DeviceType> Intrepid2::Basis< DeviceType, OutputType, PointType >::dofCoeffs_
protected

Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space spanned by the basis, := P(dofCoords_(i)) dofCoeffs_(i) are the nodal coefficients associated to basis functions i.

Rank-1 array for scalar basis with dimension (cardinality) Rank-2 array for vector basis with dimensions (cardinality, cell dimension)

Definition at line 328 of file Intrepid2_Basis.hpp.

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
OrdinalTypeArray2DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::fieldOrdinalH1PolynomialDegree_
protected

H^1 polynomial degree for each degree of freedom. Only defined for hierarchical bases right now. The number of entries per degree of freedom in this table depends on the basis type. For hypercubes, this will be the spatial dimension. We have not yet determined what this will be for simplices beyond 1D; there are not yet hierarchical simplicial bases beyond 1D in Intrepid2.

The H^1 polynomial degree is identical to the polynomial degree for H(grad) bases. For H(vol) bases, it is one higher than the polynomial degree. Since H(div) and H(curl) bases are constructed as products of H(vol) and H(grad) bases, the H^1 degree in a given dimension is the H^1 degree for the multiplicand in that dimension.

Rank-2 array with dimensions (cardinality, cell dimension)

Definition at line 350 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getFieldOrdinalsForH1Degree(), and Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getH1PolynomialDegreeOfField().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
OrdinalTypeArray2DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::fieldOrdinalPolynomialDegree_
protected

Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now. The number of entries per degree of freedom in this table depends on the basis type. For hypercubes, this will be the spatial dimension. We have not yet determined what this will be for simplices beyond 1D; there are not yet hierarchical simplicial bases beyond 1D in Intrepid2.

Rank-2 array with dimensions (cardinality, cell dimension)

Definition at line 337 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getFieldOrdinalsForDegree(), Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getPolynomialDegreeLength(), and Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getPolynomialDegreeOfField().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
OrdinalTypeArray2DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::ordinalToTag_
protected

"true" if tagToOrdinal_ and ordinalToTag_ have been initialized

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 237 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getAllDofTags(), and Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getDofTag().

template<typename DeviceType = void, typename OutputType = double, typename PointType = double>
OrdinalTypeArray3DHost Intrepid2::Basis< DeviceType, OutputType, PointType >::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 250 of file Intrepid2_Basis.hpp.

Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getAllDofOrdinal(), Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getDofCount(), and Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getDofOrdinal().


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