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

Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line. More...

#include <Intrepid2_DerivedBasis_HVOL_QUAD.hpp>

Inheritance diagram for Intrepid2::Basis_Derived_HVOL_QUAD< HVOL_LINE >:
Intrepid2::Basis_TensorBasis< HVOL_LINE::BasisBase >

Public Types

using ExecutionSpace = typename HVOL_LINE::ExecutionSpace
 
using OutputValueType = typename HVOL_LINE::OutputValueType
 
using PointValueType = typename HVOL_LINE::PointValueType
 
using OutputViewType = typename HVOL_LINE::OutputViewType
 
using PointViewType = typename HVOL_LINE::PointViewType
 
using ScalarViewType = typename HVOL_LINE::ScalarViewType
 
using BasisBase = typename HVOL_LINE::BasisBase
 
- Public Types inherited from Intrepid2::Basis_TensorBasis< HVOL_LINE::BasisBase >
using BasisBase = BasisBaseClass
 
using BasisPtr = Teuchos::RCP< HVOL_LINE::BasisBase >
 
using DeviceType = typename HVOL_LINE::BasisBase::DeviceType
 
using ExecutionSpace = typename HVOL_LINE::BasisBase::ExecutionSpace
 
using OutputValueType = typename HVOL_LINE::BasisBase::OutputValueType
 
using PointValueType = typename HVOL_LINE::BasisBase::PointValueType
 
using OrdinalTypeArray1DHost = typename HVOL_LINE::BasisBase::OrdinalTypeArray1DHost
 
using OrdinalTypeArray2DHost = typename HVOL_LINE::BasisBase::OrdinalTypeArray2DHost
 
using OutputViewType = typename HVOL_LINE::BasisBase::OutputViewType
 
using PointViewType = typename HVOL_LINE::BasisBase::PointViewType
 
using TensorBasis = Basis_TensorBasis< BasisBaseClass >
 

Public Member Functions

 Basis_Derived_HVOL_QUAD (int polyOrder_x, int polyOrder_y, const EPointType pointType=POINTTYPE_DEFAULT)
 Constructor. More...
 
 Basis_Derived_HVOL_QUAD (int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
 Constructor. More...
 
virtual const char * getName () const override
 Returns basis name. More...
 
virtual bool requireOrientation () const override
 True if orientation is required.
 
virtual OperatorTensorDecomposition getSimpleOperatorDecomposition (const EOperator &operatorType) const override
 Returns a simple decomposition of the specified operator: what operator(s) should be applied to basis1, and what operator(s) to basis2. A one-element OperatorTensorDecomposition corresponds to a single TensorData entry; a multiple-element OperatorTensorDecomposition corresponds to a VectorData object with axialComponents = false. More...
 
virtual BasisPtr< typename
Kokkos::HostSpace::device_type,
typename
BasisBase::OutputValueType,
typename
BasisBase::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...
 
- Public Member Functions inherited from Intrepid2::Basis_TensorBasis< HVOL_LINE::BasisBase >
 Basis_TensorBasis (BasisPtr basis1, BasisPtr basis2, EFunctionSpace functionSpace=FUNCTION_SPACE_MAX, const bool useShardsCellTopologyAndTags=false)
 Constructor. More...
 
void setShardsTopologyAndTags ()
 
virtual int getNumTensorialExtrusions () const override
 
ordinal_type getTensorDkEnumeration (ordinal_type dkEnum1, ordinal_type operatorOrder1, ordinal_type dkEnum2, ordinal_type operatorOrder2) const
 Given "Dk" enumeration indices for the component bases, returns a Dk enumeration index for the composite basis. More...
 
virtual OperatorTensorDecomposition getOperatorDecomposition (const EOperator operatorType) const
 Returns a full decomposition of the specified operator. (Full meaning that all TensorBasis components are expanded into their non-TensorBasis components.)
 
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...
 
void getComponentPoints (const PointViewType inputPoints, const bool attemptTensorDecomposition, PointViewType &inputPoints1, PointViewType &inputPoints2, bool &tensorDecompositionSucceeded) const
 Method to extract component points from composite points. More...
 
virtual void getDofCoords (typename HVOL_LINE::BasisBase::ScalarViewType dofCoords) const override
 Fills in spatial locations (coordinates) of degrees of freedom (nodes) on the reference cell. More...
 
virtual void getDofCoeffs (typename HVOL_LINE::BasisBase::ScalarViewType dofCoeffs) const override
 Fills in coefficients of degrees of freedom on the reference cell. More...
 
std::vector< BasisPtr > getTensorBasisComponents () const
 
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...
 
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 void getValues (OutputViewType outputValues, const EOperator operatorType, const PointViewType inputPoints1, const PointViewType inputPoints2, bool tensorPoints) const
 Evaluation of a tensor FEM basis on a reference cell; subclasses should override this. More...
 
void getValues (OutputViewType outputValues, const PointViewType inputPoints1, const EOperator operatorType1, const PointViewType inputPoints2, const EOperator operatorType2, bool tensorPoints, double weight=1.0) const
 Evaluation of a tensor FEM basis on a reference cell. More...
 

Protected Types

using LineBasis = HVOL_LINE
 
using TensorBasis = Basis_TensorBasis< typename HVOL_LINE::BasisBase >
 

Protected Attributes

std::string name_
 
ordinal_type polyOrder_x_
 
ordinal_type polyOrder_y_
 
EPointType pointType_
 
- Protected Attributes inherited from Intrepid2::Basis_TensorBasis< HVOL_LINE::BasisBase >
BasisPtr basis1_
 
BasisPtr basis2_
 
std::vector< BasisPtr > tensorComponents_
 
std::string name_
 
int numTensorialExtrusions_
 

Detailed Description

template<class HVOL_LINE>
class Intrepid2::Basis_Derived_HVOL_QUAD< HVOL_LINE >

Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line.

Definition at line 65 of file Intrepid2_DerivedBasis_HVOL_QUAD.hpp.

Constructor & Destructor Documentation

template<class HVOL_LINE >
Intrepid2::Basis_Derived_HVOL_QUAD< HVOL_LINE >::Basis_Derived_HVOL_QUAD ( int  polyOrder_x,
int  polyOrder_y,
const EPointType  pointType = POINTTYPE_DEFAULT 
)
inline

Constructor.

Parameters
[in]polyOrder_x- the polynomial order in the x dimension.
[in]polyOrder_y- the polynomial order in the y dimension.
[in]pointType- type of lattice used for creating the DoF coordinates.

Definition at line 92 of file Intrepid2_DerivedBasis_HVOL_QUAD.hpp.

References Intrepid2::Basis_TensorBasis< BasisBase >::getName().

template<class HVOL_LINE >
Intrepid2::Basis_Derived_HVOL_QUAD< HVOL_LINE >::Basis_Derived_HVOL_QUAD ( int  polyOrder,
const EPointType  pointType = POINTTYPE_DEFAULT 
)
inline

Constructor.

Parameters
[in]polyOrder- the polynomial order to use in both dimensions.
[in]pointType- type of lattice used for creating the DoF coordinates.

Definition at line 113 of file Intrepid2_DerivedBasis_HVOL_QUAD.hpp.

Member Function Documentation

template<class HVOL_LINE >
virtual BasisPtr<typename Kokkos::HostSpace::device_type, typename BasisBase::OutputValueType, typename BasisBase::PointValueType> Intrepid2::Basis_Derived_HVOL_QUAD< HVOL_LINE >::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.

Reimplemented from Intrepid2::Basis_TensorBasis< HVOL_LINE::BasisBase >.

Definition at line 152 of file Intrepid2_DerivedBasis_HVOL_QUAD.hpp.

template<class HVOL_LINE >
virtual const char* Intrepid2::Basis_Derived_HVOL_QUAD< HVOL_LINE >::getName ( ) const
inlineoverridevirtual

Returns basis name.

Returns
the name of the basis

Reimplemented from Intrepid2::Basis_TensorBasis< HVOL_LINE::BasisBase >.

Definition at line 121 of file Intrepid2_DerivedBasis_HVOL_QUAD.hpp.

template<class HVOL_LINE >
virtual OperatorTensorDecomposition Intrepid2::Basis_Derived_HVOL_QUAD< HVOL_LINE >::getSimpleOperatorDecomposition ( const EOperator &  operatorType) const
inlineoverridevirtual

Returns a simple decomposition of the specified operator: what operator(s) should be applied to basis1, and what operator(s) to basis2. A one-element OperatorTensorDecomposition corresponds to a single TensorData entry; a multiple-element OperatorTensorDecomposition corresponds to a VectorData object with axialComponents = false.

Subclasses must override this method.

Reimplemented from Intrepid2::Basis_TensorBasis< HVOL_LINE::BasisBase >.

Definition at line 131 of file Intrepid2_DerivedBasis_HVOL_QUAD.hpp.


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