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

A basis that is the direct sum of two other bases. More...

#include <Intrepid2_DirectSumBasis.hpp>

Inheritance diagram for Intrepid2::Basis_DirectSumBasis< 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 OrdinalTypeArray1DHost = typename BasisBase::OrdinalTypeArray1DHost
 
using OrdinalTypeArray2DHost = typename BasisBase::OrdinalTypeArray2DHost
 
using OutputViewType = typename BasisBase::OutputViewType
 
using PointViewType = typename BasisBase::PointViewType
 
using ScalarViewType = typename BasisBase::ScalarViewType
 

Public Member Functions

 Basis_DirectSumBasis (BasisPtr basis1, BasisPtr basis2)
 Constructor. More...
 
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 void getDofCoords (ScalarViewType dofCoords) const override
 Fills in spatial locations (coordinates) of degrees of freedom (nodes) on the reference cell. More...
 
virtual void getDofCoeffs (ScalarViewType dofCoeffs) const override
 Fills in coefficients of degrees of freedom for Lagrangian basis on the reference cell. 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 int getNumTensorialExtrusions () const override
 

Protected Attributes

BasisPtr basis1_
 
BasisPtr basis2_
 
std::string name_
 

Detailed Description

template<typename BasisBaseClass>
class Intrepid2::Basis_DirectSumBasis< BasisBaseClass >

A basis that is the direct sum of two other bases.

The direct-sum basis is ordered such that the Basis1 members come first (and in the same order as they exist in Basis1), followed by the members of Basis2, in the same order as they exist in Basis2.

The two bases must agree in their BasisType (the return value of getBasisType()).

Definition at line 67 of file Intrepid2_DirectSumBasis.hpp.

Constructor & Destructor Documentation

template<typename BasisBaseClass>
Intrepid2::Basis_DirectSumBasis< BasisBaseClass >::Basis_DirectSumBasis ( BasisPtr  basis1,
BasisPtr  basis2 
)
inline

Constructor.

Parameters
[in]basis1- the instance of Basis1
[in]basis2- the instance of Basis2

Definition at line 93 of file Intrepid2_DirectSumBasis.hpp.

Member Function Documentation

template<typename BasisBaseClass>
virtual BasisValues<OutputValueType,DeviceType> Intrepid2::Basis_DirectSumBasis< 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 default implementation employs a trivial tensor-product structure, for compatibility across all bases. Subclasses that have tensor-product structure should override. Note that only the basic exact-sequence operators are supported at the moment: VALUE, GRAD, DIV, CURL.

Definition at line 216 of file Intrepid2_DirectSumBasis.hpp.

template<typename BasisBaseClass>
virtual void Intrepid2::Basis_DirectSumBasis< BasisBaseClass >::getDofCoeffs ( ScalarViewType  dofCoeffs) const
inlineoverridevirtual

Fills in coefficients of degrees of freedom for Lagrangian basis on the reference cell.

Parameters
[out]dofCoeffs- the container into which to place the degrees of freedom.

dofCoeffs have shape (F,D) or (F) if D=1, field dimension matches the cardinality of the basis, and D is the basis dimension.

Degrees of freedom coefficients are such that (dofCoords_(j)) dofCoeffs_(j) = , where are the basis and the Kronecker delta. Note that getDofCoeffs() is supported only for Lagrangian bases.

Definition at line 305 of file Intrepid2_DirectSumBasis.hpp.

template<typename BasisBaseClass>
virtual void Intrepid2::Basis_DirectSumBasis< BasisBaseClass >::getDofCoords ( ScalarViewType  dofCoords) const
inlineoverridevirtual

Fills in spatial locations (coordinates) of degrees of freedom (nodes) on the reference cell.

Parameters
[out]dofCoords- the container into which to place the degrees of freedom.

dofCoords should have shape (F,D), where the field dimension matches the cardinality of the basis, and D is the spatial dimension of the topology on which the basis is defined.

Note that getDofCoords() is not supported by all bases; in particular, hierarchical bases do not generally support this.

Definition at line 282 of file Intrepid2_DirectSumBasis.hpp.

template<typename BasisBaseClass>
virtual const char* Intrepid2::Basis_DirectSumBasis< BasisBaseClass >::getName ( ) const
inlineoverridevirtual
template<typename BasisBaseClass>
virtual void Intrepid2::Basis_DirectSumBasis< 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 TensorBasis and DirectSumBasis and their subclasses. It allows a reduced memory footprint and optimized integration, etc.

Definition at line 346 of file Intrepid2_DirectSumBasis.hpp.

template<typename BasisBaseClass>
virtual void Intrepid2::Basis_DirectSumBasis< 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 380 of file Intrepid2_DirectSumBasis.hpp.


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