Intrepid2
Public Types | Public Member Functions | Protected Attributes | List of all members
Intrepid2::SerendipityBasis< BasisBaseClass > Class Template Reference

Serendipity Basis, defined as the sub-basis of a provided basis, consisting of basis elements for which tensorial component polynomial orders satisfy the Serendipity criterion. More...

#include <Intrepid2_SerendipityBasis.hpp>

Inheritance diagram for Intrepid2::SerendipityBasis< BasisBaseClass >:

Public Types

using BasisBase = BasisBaseClass
 
using BasisPtr = Teuchos::RCP< BasisBase >
 
using DeviceType = typename BasisBase::DeviceType
 
using ExecutionSpace = typename BasisBase::ExecutionSpace
 
using OutputValueType = typename BasisBase::OutputValueType
 
using PointValueType = typename BasisBase::PointValueType
 
using OrdinalTypeArray1D = typename BasisBase::OrdinalTypeArray1D
 
using OrdinalTypeArray1DHost = typename BasisBase::OrdinalTypeArray1DHost
 
using OrdinalTypeArray2DHost = typename BasisBase::OrdinalTypeArray2DHost
 
using OutputViewType = typename BasisBase::OutputViewType
 
using PointViewType = typename BasisBase::PointViewType
 

Public Member Functions

 SerendipityBasis (BasisPtr fullBasis)
 Constructor. More...
 
virtual int getNumTensorialExtrusions () const override
 
virtual BasisValues
< OutputValueType, DeviceType > 
allocateBasisValues (TensorPoints< PointValueType, DeviceType > points, const EOperator operatorType=OPERATOR_VALUE) const override
 Allocate BasisValues container suitable for passing to the getValues() variant that takes a TensorPoints container as argument. More...
 
virtual const char * getName () const override
 Returns basis name. More...
 
virtual void getValues (BasisValues< OutputValueType, DeviceType > outputValues, const TensorPoints< PointValueType, DeviceType > inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
 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 outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
 Evaluation of a FEM basis on a reference cell. More...
 
virtual HostBasisPtr
< 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...
 
BasisPtr getUnderlyingBasis () const
 Returns a pointer to the underlying full basis. More...
 
OrdinalTypeArray1D ordinalMap () const
 Returns the ordinal map from the Serendipity basis ordinal to the ordinal in the underlying full basis. More...
 

Protected Attributes

BasisPtr fullBasis_
 
std::string name_
 
int numTensorialExtrusions_
 
OrdinalTypeArray1D ordinalMap_
 

Detailed Description

template<typename BasisBaseClass = void>
class Intrepid2::SerendipityBasis< BasisBaseClass >

Serendipity Basis, defined as the sub-basis of a provided basis, consisting of basis elements for which tensorial component polynomial orders satisfy the Serendipity criterion.

Definition at line 20 of file Intrepid2_SerendipityBasis.hpp.

Constructor & Destructor Documentation

template<typename BasisBaseClass = void>
Intrepid2::SerendipityBasis< BasisBaseClass >::SerendipityBasis ( BasisPtr  fullBasis)
inline

Constructor.

Parameters
[in]basis- the full, hierarchical basis

Definition at line 51 of file Intrepid2_SerendipityBasis.hpp.

Member Function Documentation

template<typename BasisBaseClass = void>
virtual BasisValues<OutputValueType,DeviceType> Intrepid2::SerendipityBasis< BasisBaseClass >::allocateBasisValues ( TensorPoints< PointValueType, DeviceType >  points,
const EOperator  operatorType = OPERATOR_VALUE 
) const
inlineoverridevirtual

Allocate BasisValues container suitable for passing to the getValues() variant that takes a TensorPoints container as argument.

The basic exact-sequence operators are supported (VALUE, GRAD, DIV, CURL), as are the Dn operators (OPERATOR_D1 through OPERATOR_D10).

Definition at line 189 of file Intrepid2_SerendipityBasis.hpp.

template<typename BasisBaseClass = void>
virtual HostBasisPtr<OutputValueType, PointValueType> Intrepid2::SerendipityBasis< BasisBaseClass >::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.

Definition at line 263 of file Intrepid2_SerendipityBasis.hpp.

template<typename BasisBaseClass = void>
virtual const char* Intrepid2::SerendipityBasis< BasisBaseClass >::getName ( ) const
inlineoverridevirtual

Returns basis name.

Returns
the name of the basis

Definition at line 207 of file Intrepid2_SerendipityBasis.hpp.

template<typename BasisBaseClass = void>
BasisPtr Intrepid2::SerendipityBasis< BasisBaseClass >::getUnderlyingBasis ( ) const
inline

Returns a pointer to the underlying full basis.

Returns
Pointer to the underlying full basis.

Definition at line 274 of file Intrepid2_SerendipityBasis.hpp.

template<typename BasisBaseClass = void>
virtual void Intrepid2::SerendipityBasis< BasisBaseClass >::getValues ( BasisValues< OutputValueType, DeviceType >  outputValues,
const TensorPoints< PointValueType, DeviceType >  inputPoints,
const EOperator  operatorType = OPERATOR_VALUE 
) const
inlineoverridevirtual

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

This is the preferred getValues() method for SerendipityBasis. It allows a reduced memory footprint and optimized integration, etc.

Definition at line 224 of file Intrepid2_SerendipityBasis.hpp.

template<typename BasisBaseClass = void>
virtual void Intrepid2::SerendipityBasis< BasisBaseClass >::getValues ( OutputViewType  outputValues,
const PointViewType  inputPoints,
const EOperator  operatorType = OPERATOR_VALUE 
) const
inlineoverridevirtual

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.

Definition at line 251 of file Intrepid2_SerendipityBasis.hpp.

template<typename BasisBaseClass = void>
OrdinalTypeArray1D Intrepid2::SerendipityBasis< BasisBaseClass >::ordinalMap ( ) const
inline

Returns the ordinal map from the Serendipity basis ordinal to the ordinal in the underlying full basis.

Returns
Ordinal map from the Serendipity basis ordinal to the ordinal in the underlying full basis. (Indices to the container are Serendipity basis ordinals; values are full basis ordinals.)

Definition at line 283 of file Intrepid2_SerendipityBasis.hpp.


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