Intrepid2
|
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 (OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const |
Evaluation of a FEM basis on a reference cell. More... | |
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 °rees) 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 > | getFieldOrdinalsForDegree (std::vector< int > °rees) 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 | 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... | |
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... | |
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... | |
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... | |
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. | |
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... | |
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 Degree of freedom ordinals and tags for details.
Note the use of Kokkos::LayoutStride as the layout for various View definitions (including the argument for 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 getMatchingViewWithLabel() provided in Intrepid_Utils.hpp.
Definition at line 71 of file Intrepid2_Basis.hpp.
|
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 tensor-product structure should override. Note that only the basic exact-sequence operators are supported at the moment: VALUE, GRAD, DIV, CURL.
Definition at line 372 of file Intrepid2_Basis.hpp.
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 869 of file Intrepid2_BasisDef.hpp.
|
inline |
Retrieves all DoF tags.
Definition at line 803 of file Intrepid2_Basis.hpp.
|
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.
Definition at line 685 of file Intrepid2_Basis.hpp.
|
inline |
|
inline |
Returns cardinality of the basis.
Definition at line 656 of file Intrepid2_Basis.hpp.
|
inline |
Returns the type of coordinate system for which the basis is defined.
Definition at line 705 of file Intrepid2_Basis.hpp.
|
inline |
Returns the degree of the basis.
Definition at line 666 of file Intrepid2_Basis.hpp.
|
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)
Definition at line 515 of file Intrepid2_Basis.hpp.
|
inline |
DoF count for specified subcell.
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 |
Definition at line 717 of file Intrepid2_Basis.hpp.
|
inline |
DoF tag to ordinal lookup.
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. |
Definition at line 743 of file Intrepid2_Basis.hpp.
|
inline |
DoF ordinal to DoF tag lookup.
dofOrd | [in] - ordinal of the DoF whose tag is being retrieved |
Definition at line 785 of file Intrepid2_Basis.hpp.
Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getDofCount().
|
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.
degrees | [in] - 1D host ordinal array of length specified by getPolynomialDegreeLength(), indicating what the maximum degree in each dimension should be |
Definition at line 527 of file Intrepid2_Basis.hpp.
|
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.
degrees | [in] - std::vector<int> of length specified by getPolynomialDegreeLength(), indicating what the maximum degree in each dimension should be |
Definition at line 563 of file Intrepid2_Basis.hpp.
|
inline |
Returns the function space for the basis.
Definition at line 675 of file Intrepid2_Basis.hpp.
|
inlinevirtual |
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_type, but is otherwise identical to this.
Reimplemented in Intrepid2::IntegratedLegendreBasis_HGRAD_TET< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::IntegratedLegendreBasis_HGRAD_TRI< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::IntegratedLegendreBasis_HGRAD_LINE< DeviceType, OutputScalar, PointScalar, defineVertexFunctions, useMinusOneToOneReferenceElement >, Intrepid2::Basis_HGRAD_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::LegendreBasis_HVOL_LINE< DeviceType, OutputScalar, PointScalar >, Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_Cn_FEM< 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_HDIV_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, and Intrepid2::Basis_HVOL_C0_FEM< DeviceType, outputValueType, pointValueType >.
Definition at line 832 of file Intrepid2_Basis.hpp.
|
inlinevirtual |
Returns basis name.
Reimplemented in Intrepid2::IntegratedLegendreBasis_HGRAD_TET< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::IntegratedLegendreBasis_HGRAD_TRI< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::IntegratedLegendreBasis_HGRAD_LINE< DeviceType, OutputScalar, PointScalar, defineVertexFunctions, useMinusOneToOneReferenceElement >, Intrepid2::Basis_HGRAD_HEX_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_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_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::LegendreBasis_HVOL_LINE< DeviceType, OutputScalar, PointScalar >, Intrepid2::Basis_HVOL_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_WEDGE_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_COMP12_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_LINE_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TET_I1_FEM< 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_HGRAD_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_HEX_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_QUAD_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C2_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HVOL_QUAD_Cn_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_WEDGE_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_QUAD_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_LINE_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TRI_C1_FEM< DeviceType, outputValueType, pointValueType >, and Intrepid2::Basis_HVOL_C0_FEM< DeviceType, outputValueType, pointValueType >.
Definition at line 639 of file Intrepid2_Basis.hpp.
|
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.
fieldOrdinal | [in] - ordinal of the basis function whose polynomial degree is requested. |
Definition at line 587 of file Intrepid2_Basis.hpp.
|
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.
fieldOrdinal | [in] - ordinal of the basis function whose polynomial degree is requested. |
Definition at line 610 of file Intrepid2_Basis.hpp.
|
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).
[in] | subCellDim | - dimension of subCell |
[in] | subCellOrd | - position of the subCell among of the subCells having the same dimension |
Reimplemented in Intrepid2::IntegratedLegendreBasis_HGRAD_TET< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, Intrepid2::IntegratedLegendreBasis_HGRAD_TRI< DeviceType, OutputScalar, PointScalar, defineVertexFunctions >, 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_HDIV_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_TET_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_HEX_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HCURL_WEDGE_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_HEX_Cn_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_HCURL_TET_I1_FEM< 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_TRI_I1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HDIV_HEX_In_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_QUAD_C1_FEM< DeviceType, outputValueType, pointValueType >, Intrepid2::Basis_HGRAD_TET_C1_FEM< DeviceType, outputValueType, pointValueType >, and Intrepid2::Basis_HGRAD_TRI_C1_FEM< DeviceType, outputValueType, pointValueType >.
Definition at line 822 of file Intrepid2_Basis.hpp.
|
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.
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 |
Definition at line 422 of file Intrepid2_Basis.hpp.
|
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.
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 442 of file Intrepid2_Basis.hpp.
|
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.
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 |
Definition at line 486 of file Intrepid2_Basis.hpp.
|
inlineprotected |
Fills ordinalToTag_ and tagToOrdinal_ by basis-specific tag data.
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 263 of file Intrepid2_Basis.hpp.
|
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 325 of file Intrepid2_Basis.hpp.
|
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 334 of file Intrepid2_Basis.hpp.
|
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.
Definition at line 234 of file Intrepid2_Basis.hpp.
Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getAllDofTags().
|
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.
Definition at line 247 of file Intrepid2_Basis.hpp.
Referenced by Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getAllDofOrdinal(), and Intrepid2::Basis< DeviceType, OutputScalar, PointScalar >::getDofCount().