Intrepid2
Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI< DeviceType, outputValueType, pointValueType > Class Template Reference

Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,1] reference line cell, using Jacobi polynomials. More...

#include <Intrepid2_HGRAD_LINE_Cn_FEM_JACOBI.hpp>

Inheritance diagram for Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI< DeviceType, outputValueType, pointValueType >:
Intrepid2::Basis< void, double, double >

Public Types

typedef double value_type
using BasisBase = Basis<DeviceType,outputValueType,pointValueType>
Public Types inherited from Intrepid2::Basis< void, double, double >
using 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
 (Kokkos) Execution space for basis.
using OutputValueType
 Output value type for basis; default is double.
using PointValueType
 Point value type for basis; default is double.
using OrdinalViewType
 View type for ordinal.
using EBasisViewType
 View for basis type.
using ECoordinatesViewType
 View for coordinate system type.
using OrdinalTypeArray1DHost
 View type for 1d host array.
using OrdinalTypeArray2DHost
 View type for 2d host array.
using OrdinalTypeArray3DHost
 View type for 3d host array.
using OrdinalTypeArrayStride1DHost
 View type for 1d host array.
using OrdinalTypeArray1D
 View type for 1d device array.
using OrdinalTypeArray2D
 View type for 2d device array.
using OrdinalTypeArray3D
 View type for 3d device array.
using OrdinalTypeArrayStride1D
 View type for 1d device array.
typedef ScalarTraits< double >::scalar_type scalarType
 Scalar type for point values.
using OutputViewType
 View type for basis value output.
using PointViewType
 View type for input points.
using ScalarViewType
 View type for scalars.

Public Member Functions

 Basis_HGRAD_LINE_Cn_FEM_JACOBI (const ordinal_type order, const double alpha=0, const double beta=0)
 Constructor.
virtual void getValues (const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
virtual KOKKOS_INLINE_FUNCTION void getValues (OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
 Team-level evaluation of basis functions on a reference cell.
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.
virtual void getValues (OutputViewType, const PointViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
 Evaluation of an FVD basis evaluation on a physical cell.
Public Member Functions inherited from Intrepid2::Basis< void, double, double >
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).
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.
virtual void getScratchSpaceSize (ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
 Return the size of the scratch space, in bytes, needed for using the team-level implementation of getValues.
virtual KOKKOS_INLINE_FUNCTION void getValues (OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
 Team-level evaluation of basis functions on a reference cell.
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, \alpha_i := P(dofCoords(i)) \cdot dofCoeffs(i) are the nodal coefficients associated to basis function i.
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.
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.
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.
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.
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.
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.
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.
virtual bool requireOrientation () const
 True if orientation is required.
ordinal_type getCardinality () const
 Returns cardinality of the basis.
ordinal_type getDegree () const
 Returns the degree of the basis.
EFunctionSpace getFunctionSpace () const
 Returns the function space for the basis.
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.
EBasis getBasisType () const
 Returns the basis type.
ECoordinates getCoordinateSystem () const
 Returns the type of coordinate system for which the basis is defined.
ordinal_type getDofCount (const ordinal_type subcDim, const ordinal_type subcOrd) const
 DoF count for specified subcell.
ordinal_type getDofOrdinal (const ordinal_type subcDim, const ordinal_type subcOrd, const ordinal_type subcDofOrd) const
 DoF tag to ordinal lookup.
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.
const OrdinalTypeArray2DHost getAllDofTags () const
 Retrieves all DoF tags.
virtual BasisPtr< DeviceType, OutputValueType, PointValueTypegetSubCellRefBasis (const ordinal_type subCellDim, const ordinal_type subCellOrd) const
 returns the basis associated to a subCell.
ordinal_type getDomainDimension () const
 Returns the spatial dimension of the domain of the basis; this is equal to getBaseCellTopology().getDimension() + getNumTensorialExtrusions().
virtual HostBasisPtr< OutputValueType, PointValueTypegetHostBasis () const
 Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_type, but is otherwise identical to this.

Private Attributes

double alpha_
double beta_

Additional Inherited Members

Protected Member Functions inherited from Intrepid2::Basis< void, double, double >
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.
Protected Attributes inherited from Intrepid2::Basis< void, double, double >
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.
unsigned basisCellTopologyKey_
 Identifier of the 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 base topology, 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_
 The function space in which the basis is defined.
OrdinalTypeArray2DHost ordinalToTag_
 "true" if tagToOrdinal_ and ordinalToTag_ have been initialized
OrdinalTypeArray3DHost tagToOrdinal_
 DoF tag to ordinal lookup table.
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, \alpha_i := P(dofCoords_(i)) \cdot dofCoeffs_(i) are the nodal coefficients associated to basis functions i.
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.
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.

Detailed Description

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

Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,1] reference line cell, using Jacobi polynomials.

Implements Jacobi basis of variable order $n$ on the reference [-1,1] line cell. Jacobi polynomials depend on three parameters $ \alpha $, $ \beta $, and $ n $ and are defined via the so-called Gamma function by

\‍[P_n^{(\alpha,\beta)} (z) = 
\frac{\Gamma (\alpha+n+1)}{n!\Gamma (\alpha+\beta+n+1)}
\sum_{m=0}^n {n\choose m}
\frac{\Gamma (\alpha + \beta + n + m + 1)}{\Gamma (\alpha + m + 1)} \left(\frac{z-1}{2}\right)^m
\‍]

The basis has cardinality $n+1$ and spans a COMPLETE linear polynomial space. Basis functions are dual to a unisolvent set of degrees of freedom (DoF) enumerated as follows:

Basis order DoF tag table DoF definition
subc dim subc ordinal subc DoF tag subc num DoFs
0 1 0 0 1 $ P_0^{(\alpha,\beta)} $
1 1 0 0-1 2 $ P_0^{(\alpha,\beta)}, P_1^{(\alpha,\beta)} $
2 1 0 0-2 3 $ P_0^{(\alpha,\beta)}, P_1^{(\alpha,\beta)}, P_2^{(\alpha,\beta)} $
3 1 0 0-3 4 $ P_0^{(\alpha,\beta)}, P_1^{(\alpha,\beta)}, ..., P_3^{(\alpha,\beta)} $
... 1 0 ... ... ...
n 1 0 0-n n+1 $ P_0^{(\alpha,\beta)}, P_1^{(\alpha,\beta)}, ..., P_n^{(\alpha,\beta)} $

For example, for Legendre polynomials ( $\alpha=\beta=0$), the first 11 bases are given by

Basis order DoF tag table DoF definition
subc dim subc ordinal subc DoF tag subc num DoFs
0 1 0 0 1 $ 1 $
1 1 0 0-1 2 and: $ x $
2 1 0 0-2 3 and: $ \frac{1}{2} (3x^2-1) $
3 1 0 0-3 4 and: $ \frac{1}{2} (5x^3-3x) $
4 1 0 0-4 5 and: $ \frac{1}{8} (35x^4-30x^2+3) $
5 1 0 0-5 6 and: $ \frac{1}{8} (63x^5-70x^3+15x) $
6 1 0 0-6 7 and: $ \frac{1}{16} (231x^6-315x^4+105x^2-5) $
7 1 0 0-7 8 and: $ \frac{1}{16} (429x^7-693x^5+315x^3-35x) $
8 1 0 0-8 9 and: $ \frac{1}{128} (6435x^8-12012x^6+6930x^4-1260x^2+35) $
9 1 0 0-9 10 and: $ \frac{1}{128} (12155x^9-25740x^7+18018x^5-4620x^3+315x) $
10 1 0 0-1011 and: $ \frac{1}{128} (46189x^{10}-109395x^8+90090x^6-30030x^4+3465x^2-63) $

Definition at line 199 of file Intrepid2_HGRAD_LINE_Cn_FEM_JACOBI.hpp.

Member Typedef Documentation

◆ BasisBase

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
using Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI< DeviceType, outputValueType, pointValueType >::BasisBase = Basis<DeviceType,outputValueType,pointValueType>

Definition at line 204 of file Intrepid2_HGRAD_LINE_Cn_FEM_JACOBI.hpp.

◆ value_type

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
typedef double Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI< DeviceType, outputValueType, pointValueType >::value_type

Definition at line 202 of file Intrepid2_HGRAD_LINE_Cn_FEM_JACOBI.hpp.

Constructor & Destructor Documentation

◆ Basis_HGRAD_LINE_Cn_FEM_JACOBI()

Member Function Documentation

◆ getValues() [1/5]

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
virtual void Intrepid2::Basis< DeviceType, outputValueType, pointValueType >::getValues ( BasisValues< OutputValueType, DeviceType > outputValues,
const TensorPoints< PointValueType, DeviceType > inputPoints,
const EOperator operatorType = OPERATOR_VALUE ) const
inline

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

◆ getValues() [2/5]

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
virtual void Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI< DeviceType, outputValueType, pointValueType >::getValues ( const ExecutionSpace & space,
OutputViewType outputValues,
const PointViewType inputPoints,
const EOperator operatorType = OPERATOR_VALUE ) const
inlineoverridevirtual

Definition at line 225 of file Intrepid2_HGRAD_LINE_Cn_FEM_JACOBI.hpp.

◆ getValues() [3/5]

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
virtual void Intrepid2::Basis< DeviceType, outputValueType, pointValueType >::getValues ( OutputViewType outputValues,
const PointViewType inputPoints,
const EOperator operatorType = OPERATOR_VALUE ) const
inline

Definition at line 467 of file Intrepid2_Basis.hpp.

◆ getValues() [4/5]

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
virtual KOKKOS_INLINE_FUNCTION void Intrepid2::Basis< DeviceType, outputValueType, pointValueType >::getValues ( OutputViewType ,
const PointViewType ,
const EOperator ,
const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type & teamMember,
const typename ExecutionSpace::scratch_memory_space & scratchStorage,
const ordinal_type subcellDim = -1,
const ordinal_type subcellOrdinal = -1 ) const
inline

Team-level evaluation of basis functions on a reference cell.

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

The interface allow also to select basis functions associated to a particular entity. As an example, if subcellDim==1 (edges) and subcellOrdinal==0, outputValues will contain all the basis functions associated with the first edge. outputValues will contain all the cell basis functions when the default value (-1) is used for subcellDim and subcellOrdinal

Parameters
outputValues[out] - variable rank array with the basis values
inputPoints[in] - rank-2 array (P,D) with the evaluation points
operatorType[in] - the operator acting on the basis functions
teamMember[in] - team member of the Kokkos::TemaPolicy
scratchStorage[in] - scratch space to use by each team
subcellDim[in] - the dimension of the subcells, the default values of -1 returns basis functions associated to subcells of all dimensions
subcellOrdinal[in] - the ordinal of the subcell, the default values of -1 returns basis functions associated to subcells of all ordinals
Remarks
This function is supposed to be called within a TeamPolicy kernel. The size of the required scratch space is determined by the getScratchSpaceSize function.

Definition at line 426 of file Intrepid2_Basis.hpp.

◆ getValues() [5/5]

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
virtual void Intrepid2::Basis< DeviceType, outputValueType, pointValueType >::getValues ( OutputViewType ,
const PointViewType ,
const PointViewType ,
const EOperator = OPERATOR_VALUE ) const
inline

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

Member Data Documentation

◆ alpha_

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
double Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI< DeviceType, outputValueType, pointValueType >::alpha_
private

Definition at line 248 of file Intrepid2_HGRAD_LINE_Cn_FEM_JACOBI.hpp.

◆ beta_

template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
double Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI< DeviceType, outputValueType, pointValueType >::beta_
private

Definition at line 248 of file Intrepid2_HGRAD_LINE_Cn_FEM_JACOBI.hpp.


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