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

Implementation of the default HVOL-compatible FEM basis of degree n on Quadrilateral cell Implements Lagrangian basis of degree n on the reference Quadrilateral cell using a tensor product of points. The degrees of freedom are point evaluation at points in the interior of the Quadrilateral. More...

#include <Intrepid2_HVOL_QUAD_Cn_FEM.hpp>

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

Public Types

typedef Basis< ExecSpaceType,
outputValueType,
pointValueType >
::ordinal_type_array_1d_host 
ordinal_type_array_1d_host
 
typedef Basis< ExecSpaceType,
outputValueType,
pointValueType >
::ordinal_type_array_2d_host 
ordinal_type_array_2d_host
 
typedef Basis< ExecSpaceType,
outputValueType,
pointValueType >
::ordinal_type_array_3d_host 
ordinal_type_array_3d_host
 
typedef Basis< ExecSpaceType,
outputValueType,
pointValueType >
::outputViewType 
outputViewType
 
typedef Basis< ExecSpaceType,
outputValueType,
pointValueType >
::pointViewType 
pointViewType
 
typedef Basis< ExecSpaceType,
outputValueType,
pointValueType >
::scalarViewType 
scalarViewType
 
- Public Types inherited from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >
typedef Kokkos::View
< ordinal_type, ExecSpaceType > 
ordinal_view_type
 View type for ordinal.
 
typedef Kokkos::View< EBasis,
ExecSpaceType > 
ebasis_view_type
 View for basis type.
 
typedef Kokkos::View
< ECoordinates, ExecSpaceType > 
ecoordinates_view_type
 View for coordinate system type.
 
typedef Kokkos::View
< ordinal_type *,typename
ExecSpaceType::array_layout,
Kokkos::HostSpace > 
ordinal_type_array_1d_host
 View type for 1d host array.
 
typedef Kokkos::View
< ordinal_type **,typename
ExecSpaceType::array_layout,
Kokkos::HostSpace > 
ordinal_type_array_2d_host
 View type for 2d host array.
 
typedef Kokkos::View
< ordinal_type ***, typename
ExecSpaceType::array_layout,
Kokkos::HostSpace > 
ordinal_type_array_3d_host
 View type for 3d host array.
 
typedef Kokkos::View
< ordinal_type
*, Kokkos::LayoutStride,
Kokkos::HostSpace > 
ordinal_type_array_stride_1d_host
 View type for 1d host array.
 
typedef Kokkos::View
< ordinal_type *,ExecSpaceType > 
ordinal_type_array_1d
 View type for 1d device array.
 
typedef Kokkos::View
< ordinal_type
**,ExecSpaceType > 
ordinal_type_array_2d
 View type for 2d device array.
 
typedef Kokkos::View
< ordinal_type
***, ExecSpaceType > 
ordinal_type_array_3d
 View type for 3d device array.
 
typedef Kokkos::View
< ordinal_type
*, Kokkos::LayoutStride,
ExecSpaceType > 
ordinal_type_array_stride_1d
 View type for 1d device array.
 
typedef ScalarTraits
< pointValueType >
::scalar_type 
scalarType
 Scalar type for point values.
 
typedef Kokkos::DynRankView
< outputValueType,
Kokkos::LayoutStride,
ExecSpaceType > 
outputViewType
 View type for basis value output.
 
typedef Kokkos::DynRankView
< pointValueType,
Kokkos::LayoutStride,
ExecSpaceType > 
pointViewType
 View type for input points.
 
typedef Kokkos::DynRankView
< scalarType,
Kokkos::LayoutStride,
ExecSpaceType > 
scalarViewType
 View type for scalars.
 

Public Member Functions

 Basis_HVOL_QUAD_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
 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...
 
virtual bool requireOrientation () const
 True if orientation is required.
 
- 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...
 
ordinal_type getCardinality () const
 Returns cardinality of the basis. More...
 
ordinal_type getDegree () const
 Returns the degree of 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 ordinal_type_array_3d_host getAllDofOrdinal () const
 DoF tag to ordinal data structure.
 
const
ordinal_type_array_stride_1d_host 
getDofTag (const ordinal_type dofOrd) const
 DoF ordinal to DoF tag lookup. More...
 
const ordinal_type_array_2d_host getAllDofTags () const
 Retrieves all DoF tags. More...
 

Private Attributes

Kokkos::DynRankView< typename
scalarViewType::value_type,
ExecSpaceType > 
vinv_
 inverse of Generalized Vandermonde matrix (isotropic order)
 

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.
 
ordinal_type_array_2d_host ordinalToTag_
 "true" if tagToOrdinal_ and ordinalToTag_ have been initialized More...
 
ordinal_type_array_3d_host 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...
 

Detailed Description

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

Implementation of the default HVOL-compatible FEM basis of degree n on Quadrilateral cell Implements Lagrangian basis of degree n on the reference Quadrilateral cell using a tensor product of points. The degrees of freedom are point evaluation at points in the interior of the Quadrilateral.

Definition at line 165 of file Intrepid2_HVOL_QUAD_Cn_FEM.hpp.

Member Function Documentation

template<typename ExecSpaceType = void, typename outputValueType = double, typename pointValueType = double>
virtual void Intrepid2::Basis_HVOL_QUAD_Cn_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 222 of file Intrepid2_HVOL_QUAD_Cn_FEM.hpp.

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

template<typename ExecSpaceType = void, typename outputValueType = double, typename pointValueType = double>
virtual const char* Intrepid2::Basis_HVOL_QUAD_Cn_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 236 of file Intrepid2_HVOL_QUAD_Cn_FEM.hpp.

template<typename ExecSpaceType = void, typename outputValueType = double, typename pointValueType = double>
virtual void Intrepid2::Basis_HVOL_QUAD_Cn_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 185 of file Intrepid2_HVOL_QUAD_Cn_FEM.hpp.

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


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