Intrepid2
Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE > Class Template Reference
Inheritance diagram for Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >:
Intrepid2::Basis_TensorBasis< HGRAD_LINE::BasisBase >

Public Types

using ExecutionSpace = typename HGRAD_LINE::ExecutionSpace
using OutputValueType = typename HGRAD_LINE::OutputValueType
using PointValueType = typename HGRAD_LINE::PointValueType
using OutputViewType = typename HGRAD_LINE::OutputViewType
using PointViewType = typename HGRAD_LINE::PointViewType
using ScalarViewType = typename HGRAD_LINE::ScalarViewType
using LineBasis = HGRAD_LINE
using QuadBasis = Intrepid2::Basis_Derived_HGRAD_QUAD<HGRAD_LINE>
using BasisBase = typename HGRAD_LINE::BasisBase
using TensorBasis = Basis_TensorBasis<BasisBase>
Public Types inherited from Intrepid2::Basis_TensorBasis< HGRAD_LINE::BasisBase >
using BasisBase
using BasisPtr
using DeviceType
using ExecutionSpace
using OutputValueType
using PointValueType
using OrdinalTypeArray1DHost
using OrdinalTypeArray2DHost
using OutputViewType
using PointViewType
using TensorBasis

Public Member Functions

 Basis_Derived_HGRAD_HEX (int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT)
 Constructor.
 Basis_Derived_HGRAD_HEX (int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
 Constructor.
virtual bool requireOrientation () const override
 True if orientation is required.
virtual const char * getName () const override
 Returns basis name.
Teuchos::RCP< BasisBase > getSubCellRefBasis (const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
 returns the basis associated to a subCell.
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.
Public Member Functions inherited from Intrepid2::Basis_TensorBasis< HGRAD_LINE::BasisBase >
 Basis_TensorBasis (BasisPtr basis1, BasisPtr basis2, EFunctionSpace functionSpace=FUNCTION_SPACE_MAX, const bool useShardsCellTopologyAndTags=false)
 Constructor.
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.
virtual OperatorTensorDecomposition getSimpleOperatorDecomposition (const EOperator &operatorType) const
 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.
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.
void getComponentPoints (const PointViewType inputPoints, const bool attemptTensorDecomposition, PointViewType &inputPoints1, PointViewType &inputPoints2, bool &tensorDecompositionSucceeded) const
 Method to extract component points from composite points.
virtual void getDofCoords (typename BasisBase::ScalarViewType dofCoords) const override
 Fills in spatial locations (coordinates) of degrees of freedom (nodes) on the reference cell.
virtual void getDofCoeffs (typename BasisBase::ScalarViewType dofCoeffs) const override
 Fills in coefficients of degrees of freedom on the reference cell.
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.

Public Attributes

std::string name_
ordinal_type order_x_
ordinal_type order_y_
ordinal_type order_z_
EPointType pointType_

Additional Inherited Members

Protected Attributes inherited from Intrepid2::Basis_TensorBasis< HGRAD_LINE::BasisBase >
BasisPtr basis1_
BasisPtr basis2_
std::vector< BasisPtr > tensorComponents_
std::string name_
int numTensorialExtrusions_

Detailed Description

template<class HGRAD_LINE>
class Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >

Definition at line 34 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

Member Typedef Documentation

◆ BasisBase

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::BasisBase = typename HGRAD_LINE::BasisBase

Definition at line 48 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ ExecutionSpace

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::ExecutionSpace = typename HGRAD_LINE::ExecutionSpace

Definition at line 38 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ LineBasis

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::LineBasis = HGRAD_LINE

Definition at line 46 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ OutputValueType

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::OutputValueType = typename HGRAD_LINE::OutputValueType

Definition at line 39 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ OutputViewType

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::OutputViewType = typename HGRAD_LINE::OutputViewType

Definition at line 42 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ PointValueType

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::PointValueType = typename HGRAD_LINE::PointValueType

Definition at line 40 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ PointViewType

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::PointViewType = typename HGRAD_LINE::PointViewType

Definition at line 43 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ QuadBasis

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::QuadBasis = Intrepid2::Basis_Derived_HGRAD_QUAD<HGRAD_LINE>

Definition at line 47 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ ScalarViewType

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::ScalarViewType = typename HGRAD_LINE::ScalarViewType

Definition at line 44 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ TensorBasis

template<class HGRAD_LINE>
using Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::TensorBasis = Basis_TensorBasis<BasisBase>

Definition at line 49 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

Constructor & Destructor Documentation

◆ Basis_Derived_HGRAD_HEX() [1/2]

template<class HGRAD_LINE>
Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::Basis_Derived_HGRAD_HEX ( int polyOrder_x,
int polyOrder_y,
int polyOrder_z,
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]polyOrder_z- the polynomial order in the z dimension.
[in]pointType- type of lattice used for creating the DoF coordinates.

Definition at line 63 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

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

Referenced by Basis_Derived_HGRAD_HEX(), and getHostBasis().

◆ Basis_Derived_HGRAD_HEX() [2/2]

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

Constructor.

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

Definition at line 87 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

References Basis_Derived_HGRAD_HEX().

Member Function Documentation

◆ getHostBasis()

template<class HGRAD_LINE>
virtual HostBasisPtr< OutputValueType, PointValueType > Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_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< HGRAD_LINE::BasisBase >.

Definition at line 163 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

References Basis_Derived_HGRAD_HEX().

◆ getName()

template<class HGRAD_LINE>
virtual const char * Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::getName ( ) const
inlineoverridevirtual

Returns basis name.

Returns
the name of the basis

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

Definition at line 105 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ getSubCellRefBasis()

template<class HGRAD_LINE>
Teuchos::RCP< BasisBase > Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::getSubCellRefBasis ( const ordinal_type subCellDim,
const ordinal_type subCellOrd ) const
inlineoverride

returns the basis associated to a subCell.

The bases of the subCell should be the restriction to the subCell of the bases of the parent cell. TODO: test this method when different orders are used in different directions

Parameters
[in]subCellDim- dimension of subCell
[in]subCellOrd- position of the subCell among of the subCells having the same dimension
Returns
pointer to the subCell basis of dimension subCellDim and position subCellOrd

Definition at line 119 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ requireOrientation()

template<class HGRAD_LINE>
virtual bool Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::requireOrientation ( ) const
inlineoverridevirtual

True if orientation is required.

Definition at line 93 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

Member Data Documentation

◆ name_

template<class HGRAD_LINE>
std::string Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::name_

Definition at line 51 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ order_x_

template<class HGRAD_LINE>
ordinal_type Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::order_x_

Definition at line 52 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ order_y_

template<class HGRAD_LINE>
ordinal_type Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::order_y_

Definition at line 53 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ order_z_

template<class HGRAD_LINE>
ordinal_type Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::order_z_

Definition at line 54 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.

◆ pointType_

template<class HGRAD_LINE>
EPointType Intrepid2::Basis_Derived_HGRAD_HEX< HGRAD_LINE >::pointType_

Definition at line 55 of file Intrepid2_DerivedBasis_HGRAD_HEX.hpp.


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