Intrepid2
|
Implementation of the default H(div)-compatible FEM basis on Hexahedron cell. More...
#include <Intrepid2_HDIV_HEX_In_FEM.hpp>
Public Types | |
using | OrdinalTypeArray1DHost = typename Basis< ExecSpaceType, outputValueType, pointValueType >::OrdinalTypeArray1DHost |
using | OrdinalTypeArray2DHost = typename Basis< ExecSpaceType, outputValueType, pointValueType >::OrdinalTypeArray2DHost |
using | OrdinalTypeArray3DHost = typename Basis< ExecSpaceType, outputValueType, pointValueType >::OrdinalTypeArray3DHost |
using | OutputViewType = typename Basis< ExecSpaceType, outputValueType, pointValueType >::OutputViewType |
using | PointViewType = typename Basis< ExecSpaceType, outputValueType, pointValueType >::PointViewType |
using | ScalarViewType = typename Basis< ExecSpaceType, outputValueType, pointValueType >::ScalarViewType |
Public Types inherited from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType > | |
using | ExecutionSpace = ExecSpaceType |
(Kokkos) Execution space for basis. | |
using | OutputValueType = outputValueType |
Output value type for basis; default is double. | |
using | PointValueType = pointValueType |
Point value type for basis; default is double. | |
using | OrdinalViewType = Kokkos::View< ordinal_type, ExecSpaceType > |
View type for ordinal. | |
using | EBasisViewType = Kokkos::View< EBasis, ExecSpaceType > |
View for basis type. | |
using | ECoordinatesViewType = Kokkos::View< ECoordinates, ExecSpaceType > |
View for coordinate system type. | |
using | OrdinalTypeArray1DHost = Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > |
View type for 1d host array. | |
using | OrdinalTypeArray2DHost = Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > |
View type for 2d host array. | |
using | OrdinalTypeArray3DHost = Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > |
View type for 3d host array. | |
using | OrdinalTypeArrayStride1DHost = Kokkos::View< ordinal_type *, Kokkos::LayoutStride, Kokkos::HostSpace > |
View type for 1d host array. | |
using | OrdinalTypeArray1D = Kokkos::View< ordinal_type *, ExecSpaceType > |
View type for 1d device array. | |
using | OrdinalTypeArray2D = Kokkos::View< ordinal_type **, ExecSpaceType > |
View type for 2d device array. | |
using | OrdinalTypeArray3D = Kokkos::View< ordinal_type ***, ExecSpaceType > |
View type for 3d device array. | |
using | OrdinalTypeArrayStride1D = Kokkos::View< ordinal_type *, Kokkos::LayoutStride, ExecSpaceType > |
View type for 1d device array. | |
typedef ScalarTraits < pointValueType > ::scalar_type | scalarType |
Scalar type for point values. | |
using | OutputViewType = Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > |
View type for basis value output. | |
using | PointViewType = Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > |
View type for input points. | |
using | ScalarViewType = Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > |
View type for scalars. | |
Public Member Functions | |
Basis_HDIV_HEX_In_FEM (const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED) | |
Constructor. More... | |
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 |
virtual void | getDofCoeffs (ScalarViewType dofCoeffs) const |
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... | |
virtual void | getDofCoords (ScalarViewType) const |
Returns spatial locations (coordinates) of degrees of freedom on the reference cell. | |
virtual void | getDofCoeffs (ScalarViewType) 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... | |
OrdinalTypeArray1DHost | getFieldOrdinalsForDegree (OrdinalTypeArray1DHost °rees) const |
For hierarchical bases, returns the field ordinals that have at most the specified degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees. More... | |
std::vector< int > | getFieldOrdinalsForDegree (std::vector< int > °rees) const |
For hierarchical bases, returns the field ordinals that have at most the specified degree in each dimension. Assuming that these are less than or equal to the polynomial orders provided at Basis construction, the corresponding polynomials will form a superset of the Basis of the same type constructed with polynomial orders corresponding to the specified degrees. More... | |
OrdinalTypeArray1DHost | getPolynomialDegreeOfField (int fieldOrdinal) const |
For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array. More... | |
std::vector< int > | getPolynomialDegreeOfFieldAsVector (int fieldOrdinal) const |
For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spatial dimensions) for the specified basis ordinal as a host array. More... | |
int | getPolynomialDegreeLength () const |
For hierarchical bases, returns the number of entries required to specify the polynomial degree of a basis function. | |
ordinal_type | getCardinality () const |
Returns cardinality of the basis. More... | |
ordinal_type | getDegree () const |
Returns the degree of the basis. More... | |
EFunctionSpace | getFunctionSpace () const |
Returns the function space for 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 OrdinalTypeArray3DHost | getAllDofOrdinal () const |
DoF tag to ordinal data structure. | |
const OrdinalTypeArrayStride1DHost | getDofTag (const ordinal_type dofOrd) const |
DoF ordinal to DoF tag lookup. More... | |
const OrdinalTypeArray2DHost | getAllDofTags () const |
Retrieves all DoF tags. More... | |
Private Attributes | |
Kokkos::DynRankView< typename ScalarViewType::value_type, ExecSpaceType > | vinvLine_ |
inverse of Generalized Vandermonde matrix (isotropic order) | |
Kokkos::DynRankView< typename ScalarViewType::value_type, ExecSpaceType > | vinvBubble_ |
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. | |
EFunctionSpace | functionSpace_ = FUNCTION_SPACE_MAX |
The function space in which the basis is defined. | |
OrdinalTypeArray2DHost | ordinalToTag_ |
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized More... | |
OrdinalTypeArray3DHost | 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... | |
OrdinalTypeArray2DHost | fieldOrdinalPolynomialDegree_ |
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now. The number of entries per degree of freedom in this table depends on the basis type. For hypercubes, this will be the spatial dimension. We have not yet determined what this will be for simplices beyond 1D; there are not yet hierarchical simplicial bases beyond 1D in Intrepid2. More... | |
Implementation of the default H(div)-compatible FEM basis on Hexahedron cell.
Implements Raviart-Thomas basis of degree n on the reference Hexahedron cell. The basis has cardinality 3(n+1)n^2 and spans a INCOMPLETE polynomial space.
Definition at line 170 of file Intrepid2_HDIV_HEX_In_FEM.hpp.
Intrepid2::Basis_HDIV_HEX_In_FEM< SpT, OT, PT >::Basis_HDIV_HEX_In_FEM | ( | const ordinal_type | order, |
const EPointType | pointType = POINTTYPE_EQUISPACED |
||
) |
Constructor.
Product rule: z -> y -> x, x-normal first
Definition at line 324 of file Intrepid2_HDIV_HEX_In_FEMDef.hpp.
References Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getCardinality(), Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getDofTag(), and Intrepid2::Parameters::MaxOrder.
|
inlinevirtual |
Returns basis name.
Reimplemented from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >.
Definition at line 246 of file Intrepid2_HDIV_HEX_In_FEM.hpp.
|
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.
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 |
Reimplemented from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >.
Definition at line 190 of file Intrepid2_HDIV_HEX_In_FEM.hpp.
References Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getBaseCellTopology(), Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getCardinality(), and Intrepid2::Basis_HDIV_HEX_In_FEM< ExecSpaceType, outputValueType, pointValueType >::vinvLine_.