Intrepid2
Public Types | Public Member Functions | Private Attributes | List of all members
Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType > Class Template Reference

Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell. More...

#include <Intrepid2_HGRAD_TRI_Cn_FEM.hpp>

Inheritance diagram for Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >:
Intrepid2::Basis< DeviceType, outputValueType, pointValueType >

Public Types

using OrdinalTypeArray1DHost = typename Basis< DeviceType, outputValueType, pointValueType >::OrdinalTypeArray1DHost
 
using OrdinalTypeArray2DHost = typename Basis< DeviceType, outputValueType, pointValueType >::OrdinalTypeArray2DHost
 
using OrdinalTypeArray3DHost = typename Basis< DeviceType, outputValueType, pointValueType >::OrdinalTypeArray3DHost
 
using OutputViewType = typename Basis< DeviceType, outputValueType, pointValueType >::OutputViewType
 
using PointViewType = typename Basis< DeviceType, outputValueType, pointValueType >::PointViewType
 
using ScalarViewType = typename Basis< DeviceType, outputValueType, pointValueType >::ScalarViewType
 
typedef Basis< DeviceType, outputValueType, pointValueType >::scalarType scalarType
 
- Public Types inherited from Intrepid2::Basis< DeviceType, outputValueType, pointValueType >
using DeviceType = DeviceType
 (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

 Basis_HGRAD_TRI_Cn_FEM (const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
 Constructor.
 
virtual void getValues (OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
 
virtual void getDofCoords (ScalarViewType dofCoords) const override
 
virtual void getDofCoeffs (ScalarViewType dofCoeffs) const override
 
virtual const char * getName () const override
 Returns basis name. More...
 
virtual bool requireOrientation () const override
 True if orientation is required.
 
void getVandermondeInverse (ScalarViewType vinv) const
 
Kokkos::DynRankView< typename ScalarViewType::const_value_type, DeviceTypegetVandermondeInverse () const
 
ordinal_type getWorkSizePerPoint (const EOperator operatorType) const
 
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis (const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
 returns the basis associated to a subCell. More...
 
BasisPtr< typename Kokkos::HostSpace::device_type, outputValueType, pointValueType > getHostBasis () const override
 Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_type, but is otherwise identical to this. More...
 
- Public Member Functions inherited from Intrepid2::Basis< DeviceType, outputValueType, pointValueType >
OutputValueType getDummyOutputValue ()
 Dummy array to receive input arguments.
 
PointValueType getDummyPointValue ()
 Dummy array to receive input arguments.
 
Kokkos::DynRankView< OutputValueType, DeviceTypeallocateOutputView (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, DeviceTypeallocateBasisValues (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 &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 > 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...
 
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.
 
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...
 

Private Attributes

Kokkos::DynRankView< scalarType, DeviceTypevinv_
 inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the nodal basis in terms of phis_
 
EPointType pointType_
 type of lattice used for creating the DoF coordinates
 

Additional Inherited Members

- Protected Member Functions inherited from Intrepid2::Basis< DeviceType, outputValueType, pointValueType >
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 inherited from Intrepid2::Basis< DeviceType, outputValueType, pointValueType >
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_
 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, DeviceTypedofCoords_
 Coordinates of degrees-of-freedom for basis functions defined in physical space.
 
Kokkos::DynRankView< scalarType, DeviceTypedofCoeffs_
 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...
 

Detailed Description

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
class Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >

Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell.

Implements Lagrangian basis of degree n on the reference Triangle cell. The basis has cardinality (n+1)(n+2)/2 and spans a COMPLETE polynomial space of degree n. Basis functions are dual to a unisolvent set of degrees-of-freedom (DoF) defined on a lattice of order n (see PointTools). In particular, the degrees of freedom are point evaluation at

The distribution of these points is specified by the pointType argument to the class constructor. Currently, either equispaced lattice points or Warburton's warp-blend points are available.

The dof are enumerated according to the ordering on the lattice (see PointTools). In particular, dof number 0 is at the bottom left vertex (0,0). The dof increase along the lattice with points along the lines of constant x adjacent in the enumeration.

Definition at line 183 of file Intrepid2_HGRAD_TRI_Cn_FEM.hpp.

Member Function Documentation

◆ getHostBasis()

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType> Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >::getHostBasis ( ) const
inlineoverridevirtual

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 from Intrepid2::Basis< DeviceType, outputValueType, pointValueType >.

Definition at line 320 of file Intrepid2_HGRAD_TRI_Cn_FEM.hpp.

References Intrepid2::Basis< DeviceType, outputValueType, pointValueType >::basisDegree_.

◆ getName()

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
virtual const char* Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >::getName ( ) const
inlineoverridevirtual

Returns basis name.

Returns
the name of the basis

Reimplemented from Intrepid2::Basis< DeviceType, outputValueType, pointValueType >.

Definition at line 267 of file Intrepid2_HGRAD_TRI_Cn_FEM.hpp.

◆ getSubCellRefBasis()

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
BasisPtr<DeviceType,outputValueType,pointValueType> Intrepid2::Basis_HGRAD_TRI_Cn_FEM< DeviceType, outputValueType, pointValueType >::getSubCellRefBasis ( const ordinal_type  subCellDim,
const ordinal_type  subCellOrd 
) const
inlineoverridevirtual

returns the basis associated to a subCell.

The bases of the subCell are the restriction to the subCell of the bases of the parent cell.

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 from Intrepid2::Basis< DeviceType, outputValueType, pointValueType >.

Definition at line 310 of file Intrepid2_HGRAD_TRI_Cn_FEM.hpp.

References Intrepid2::Basis< DeviceType, outputValueType, pointValueType >::basisDegree_.


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