|
Intrepid2
|
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell. More...
#include <Intrepid2_HGRAD_HEX_C2_FEM.hpp>
Public Types | |
|
typedef Basis< ExecSpaceType, outputValueType, pointValueType > ::ordinal_type_array_1d_host | ordinal_type_array_1d_host |
|
typedef Basis< ExecSpaceType, outputValueType, pointValueType > ::ordinal_type_array_2d_host | ordinal_type_array_2d_host |
|
typedef Basis< ExecSpaceType, outputValueType, pointValueType > ::ordinal_type_array_3d_host | ordinal_type_array_3d_host |
|
typedef Basis< ExecSpaceType, outputValueType, pointValueType > ::outputViewType | outputViewType |
|
typedef Basis< ExecSpaceType, outputValueType, pointValueType > ::pointViewType | pointViewType |
|
typedef Basis< ExecSpaceType, outputValueType, pointValueType > ::scalarViewType | scalarViewType |
Public Types inherited from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType > | |
|
typedef Kokkos::View < ordinal_type, ExecSpaceType > | ordinal_view_type |
| View type for ordinal. | |
|
typedef Kokkos::View< EBasis, ExecSpaceType > | ebasis_view_type |
| View for basis type. | |
|
typedef Kokkos::View < ECoordinates, ExecSpaceType > | ecoordinates_view_type |
| View for coordinate system type. | |
|
typedef Kokkos::View < ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > | ordinal_type_array_1d_host |
| View type for 1d host array. | |
|
typedef Kokkos::View < ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > | ordinal_type_array_2d_host |
| View type for 2d host array. | |
|
typedef Kokkos::View < ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > | ordinal_type_array_3d_host |
| View type for 3d host array. | |
|
typedef Kokkos::View < ordinal_type *, Kokkos::LayoutStride, Kokkos::HostSpace > | ordinal_type_array_stride_1d_host |
| View type for 1d host array. | |
|
typedef Kokkos::View < ordinal_type *,ExecSpaceType > | ordinal_type_array_1d |
| View type for 1d device array. | |
|
typedef Kokkos::View < ordinal_type **,ExecSpaceType > | ordinal_type_array_2d |
| View type for 2d device array. | |
|
typedef Kokkos::View < ordinal_type ***, ExecSpaceType > | ordinal_type_array_3d |
| View type for 3d device array. | |
|
typedef Kokkos::View < ordinal_type *, Kokkos::LayoutStride, ExecSpaceType > | ordinal_type_array_stride_1d |
| View type for 1d device array. | |
|
typedef ScalarTraits < pointValueType > ::scalar_type | scalarType |
| Scalar type for point values. | |
|
typedef Kokkos::DynRankView < outputValueType, Kokkos::LayoutStride, ExecSpaceType > | outputViewType |
| View type for basis value output. | |
|
typedef Kokkos::DynRankView < pointValueType, Kokkos::LayoutStride, ExecSpaceType > | pointViewType |
| View type for input points. | |
|
typedef Kokkos::DynRankView < scalarType, Kokkos::LayoutStride, ExecSpaceType > | scalarViewType |
| View type for scalars. | |
Public Member Functions | |
| Basis_HGRAD_HEX_C2_FEM () | |
| Constructor. | |
| 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 |
| Returns spatial locations (coordinates) of degrees of freedom on the reference cell. | |
| virtual void | getDofCoeffs (scalarViewType dofCoeffs) 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... | |
| virtual const char * | getName () const |
| Returns basis name. More... | |
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 bool | requireOrientation () const |
| True if orientation is required. | |
| ordinal_type | getCardinality () const |
| Returns cardinality of the basis. More... | |
| ordinal_type | getDegree () const |
| Returns the degree of 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 ordinal_type_array_3d_host | getAllDofOrdinal () const |
| DoF tag to ordinal data structure. | |
| const ordinal_type_array_stride_1d_host | getDofTag (const ordinal_type dofOrd) const |
| DoF ordinal to DoF tag lookup. More... | |
| const ordinal_type_array_2d_host | getAllDofTags () const |
| Retrieves all DoF tags. More... | |
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. | |
| ordinal_type_array_2d_host | ordinalToTag_ |
| "true" if tagToOrdinal_ and ordinalToTag_ have been initialized More... | |
| ordinal_type_array_3d_host | 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... | |
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell.
Implements Lagrangian basis of degree 2 on the reference Hexahedron cell. The basis has cardinality 27 and spans a COMPLETE tri-quadratic polynomial space. Basis functions are dual to a unisolvent set of degrees-of-freedom (DoF) defined and enumerated as follows:
================================================================================================= | | degree-of-freedom-tag table | | | DoF |----------------------------------------------------------| DoF definition | | ordinal | subc dim | subc ordinal | subc DoF ord |subc num DoF | | |=========|==============|==============|==============|=============|===========================| | 0 | 0 | 0 | 0 | 1 | L_0(u) = u(-1,-1,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 1 | 0 | 1 | 0 | 1 | L_1(u) = u( 1,-1,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 2 | 0 | 2 | 0 | 1 | L_2(u) = u( 1, 1,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 3 | 0 | 3 | 0 | 1 | L_3(u) = u(-1, 1,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 4 | 0 | 4 | 0 | 1 | L_4(u) = u(-1,-1, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 5 | 0 | 5 | 0 | 1 | L_5(u) = u( 1,-1, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 6 | 0 | 6 | 0 | 1 | L_6(u) = u( 1, 1, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 7 | 0 | 7 | 0 | 1 | L_7(u) = u(-1, 1, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| |---------|--------------|--------------|--------------|-------------|---------------------------| | 8 | 1 | 0 | 0 | 1 | L_8(u) = u( 0,-1,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 9 | 1 | 1 | 0 | 1 | L_9(u) = u( 1, 0,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 10 | 1 | 2 | 0 | 1 | L_10(u) = u( 0, 1,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 11 | 1 | 3 | 0 | 1 | L_11(u) = u(-1, 0,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 12 | 1 | 8 | 0 | 1 | L_12(u) = u(-1,-1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 13 | 1 | 9 | 0 | 1 | L_13(u) = u( 1,-1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 14 | 1 | 10 | 0 | 1 | L_14(u) = u( 1, 1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 15 | 1 | 11 | 0 | 1 | L_15(u) = u(-1, 1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 16 | 1 | 4 | 0 | 1 | L_16(u) = u( 0,-1, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 17 | 1 | 5 | 0 | 1 | L_17(u) = u( 1, 0, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 18 | 1 | 6 | 0 | 1 | L_18(u) = u( 0, 1, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 19 | 1 | 7 | 0 | 1 | L_19(u) = u(-1, 0, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| |---------|--------------|--------------|--------------|-------------|---------------------------| | 20 | 3 | 0 | 0 | 1 | L_20(u) = u( 0, 0, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| |---------|--------------|--------------|--------------|-------------|---------------------------| | 21 | 2 | 4 | 0 | 1 | L_21(u) = u( 0, 0,-1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 22 | 2 | 5 | 0 | 1 | L_22(u) = u( 0, 0, 1) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 23 | 2 | 3 | 0 | 1 | L_23(u) = u(-1, 0, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 24 | 2 | 1 | 0 | 1 | L_24(u) = u( 1, 0, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 25 | 2 | 0 | 0 | 1 | L_25(u) = u( 0,-1, 0) | |---------|--------------|--------------|--------------|-------------|---------------------------| | 26 | 2 | 2 | 0 | 1 | L_26(u) = u( 0, 1, 0) | |=========|==============|==============|==============|=============|===========================| | MAX | maxScDim=2 | maxScOrd=12 | maxDfOrd=0 | - | | |=========|==============|==============|==============|=============|===========================|
Definition at line 222 of file Intrepid2_HGRAD_HEX_C2_FEM.hpp.
|
inlinevirtual |
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.
Rank-1 array for scalar basis with dimension (cardinality) Rank-2 array for vector basis with dimensions (cardinality, cell dimension)
Reimplemented from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >.
Definition at line 276 of file Intrepid2_HGRAD_HEX_C2_FEM.hpp.
References Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getCardinality().
|
inlinevirtual |
Returns basis name.
Reimplemented from Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >.
Definition at line 290 of file Intrepid2_HGRAD_HEX_C2_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 240 of file Intrepid2_HGRAD_HEX_C2_FEM.hpp.
References Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getBaseCellTopology(), and Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getCardinality().
1.8.5