Intrepid2
Public Types | Public Member Functions | List of all members
Intrepid2::Basis_HGRAD_TRI_C1_FEM< ExecSpaceType, outputValueType, pointValueType > Class Template Reference

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

#include <Intrepid2_HGRAD_TRI_C1_FEM.hpp>

Inheritance diagram for Intrepid2::Basis_HGRAD_TRI_C1_FEM< ExecSpaceType, outputValueType, pointValueType >:
Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >

Public Types

using OrdinalTypeArray1DHost = typename Basis< ExecSpaceType, outputValueType, pointValueType >::OrdinalTypeArray1DHost
 
using OrdinalTypeArray2DHost = typename Basis< ExecSpaceType, outputValueType, pointValueType >::OrdinalTypeArray2DHost
 
using OrdinalTypeArray3DHost = typename Basis< ExecSpaceType, outputValueType, pointValueType >::OrdinalTypeArray3DHost
 
using OutputViewType = typename Basis< ExecSpaceType, outputValueType, pointValueType >::OutputViewType
 
using PointViewType = typename Basis< ExecSpaceType, outputValueType, pointValueType >::PointViewType
 
using ScalarViewType = typename Basis< ExecSpaceType, outputValueType, pointValueType >::ScalarViewType
 
- Public Types inherited from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >
using ExecutionSpace = ExecSpaceType
 (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, ExecSpaceType >
 View type for ordinal.
 
using EBasisViewType = Kokkos::View< EBasis, ExecSpaceType >
 View for basis type.
 
using ECoordinatesViewType = Kokkos::View< ECoordinates, ExecSpaceType >
 View for coordinate system type.
 
using OrdinalTypeArray1DHost = Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace >
 View type for 1d host array.
 
using OrdinalTypeArray2DHost = Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace >
 View type for 2d host array.
 
using OrdinalTypeArray3DHost = Kokkos::View< ordinal_type ***, typename ExecSpaceType::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 *, ExecSpaceType >
 View type for 1d device array.
 
using OrdinalTypeArray2D = Kokkos::View< ordinal_type **, ExecSpaceType >
 View type for 2d device array.
 
using OrdinalTypeArray3D = Kokkos::View< ordinal_type ***, ExecSpaceType >
 View type for 3d device array.
 
using OrdinalTypeArrayStride1D = Kokkos::View< ordinal_type *, Kokkos::LayoutStride, ExecSpaceType >
 View type for 1d device array.
 
typedef ScalarTraits
< pointValueType >
::scalar_type 
scalarType
 Scalar type for point values.
 
using OutputViewType = Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType >
 View type for basis value output.
 
using PointViewType = Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType >
 View type for input points.
 
using ScalarViewType = Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType >
 View type for scalars.
 

Public Member Functions

 Basis_HGRAD_TRI_C1_FEM ()
 Constructor.
 
virtual void getValues (OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
 Evaluation of a FEM basis on a reference cell. More...
 
virtual void getDofCoords (ScalarViewType dofCoords) const
 Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
 
virtual void getDofCoeffs (ScalarViewType dofCoeffs) 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...
 
virtual const char * getName () const
 Returns basis name. More...
 
- Public Member Functions inherited from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >
OutputValueType getDummyOutputValue ()
 Dummy array to receive input arguments.
 
PointValueType getDummyPointValue ()
 Dummy array to receive input arguments.
 
virtual void getValues (OutputViewType, const PointViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
 Evaluation of an FVD basis evaluation on a physical cell. 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.
 
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...
 

Additional Inherited Members

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

Detailed Description

template<typename ExecSpaceType = void, typename outputValueType = double, typename pointValueType = double>
class Intrepid2::Basis_HGRAD_TRI_C1_FEM< ExecSpaceType, outputValueType, pointValueType >

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

Implements Lagrangian basis of degree 1 on the reference Triangle cell. The basis has cardinality 3 and spans a COMPLETE linear polynomial space. Basis functions are dual to a unisolvent set of degrees-of-freedom (DoF) defined and enumerated as follows:

=================================================================================================
|         |           degree-of-freedom-tag table                    |                           |
|   DoF   |----------------------------------------------------------|      DoF definition       |
| ordinal |  subc dim    | subc ordinal | subc DoF ord |subc num DoF |                           |
|=========|==============|==============|==============|=============|===========================|
|    0    |       0      |       0      |       0      |      1      |   L_0(u) = u(0,0)         |
|---------|--------------|--------------|--------------|-------------|---------------------------|
|    1    |       0      |       1      |       0      |      1      |   L_1(u) = u(1,0)         |
|---------|--------------|--------------|--------------|-------------|---------------------------|
|    2    |       0      |       2      |       0      |      1      |   L_2(u) = u(0,1)         |
|=========|==============|==============|==============|=============|===========================|
|   MAX   |  maxScDim=0  |  maxScOrd=2  |  maxDfOrd=0  |     -       |                           |
|=========|==============|==============|==============|=============|===========================|

Definition at line 157 of file Intrepid2_HGRAD_TRI_C1_FEM.hpp.

Member Function Documentation

template<typename ExecSpaceType = void, typename outputValueType = double, typename pointValueType = double>
virtual void Intrepid2::Basis_HGRAD_TRI_C1_FEM< ExecSpaceType, outputValueType, pointValueType >::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 from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >.

Definition at line 211 of file Intrepid2_HGRAD_TRI_C1_FEM.hpp.

References Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getCardinality().

template<typename ExecSpaceType = void, typename outputValueType = double, typename pointValueType = double>
virtual const char* Intrepid2::Basis_HGRAD_TRI_C1_FEM< ExecSpaceType, outputValueType, pointValueType >::getName ( ) const
inlinevirtual

Returns basis name.

Returns
the name of the basis

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

Definition at line 225 of file Intrepid2_HGRAD_TRI_C1_FEM.hpp.

template<typename ExecSpaceType = void, typename outputValueType = double, typename pointValueType = double>
virtual void Intrepid2::Basis_HGRAD_TRI_C1_FEM< ExecSpaceType, outputValueType, pointValueType >::getValues ( 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
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 from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >.

Definition at line 175 of file Intrepid2_HGRAD_TRI_C1_FEM.hpp.

References Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getBaseCellTopology(), and Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getCardinality().


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