|
Intrepid
|
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Wedge cell. More...
#include <Intrepid_HCURL_WEDGE_I1_FEM.hpp>
Public Member Functions | |
| Basis_HCURL_WEDGE_I1_FEM () | |
| Constructor. | |
| void | getValues (ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const |
| Evaluation of a FEM basis on a reference Wedge cell. More... | |
| void | getValues (ArrayScalar &outputValues, const ArrayScalar &inputPoints, const ArrayScalar &cellVertices, const EOperator operatorType=OPERATOR_VALUE) const |
| FVD basis evaluation: invocation of this method throws an exception. | |
Public Member Functions inherited from Intrepid::Basis< Scalar, ArrayScalar > | |
| virtual | ~Basis () |
| Destructor. | |
| virtual int | getCardinality () const |
| Returns cardinality of the basis. More... | |
| virtual int | getDegree () const |
| Returns the degree of the basis. More... | |
| virtual const shards::CellTopology | getBaseCellTopology () const |
| Returns the base cell topology for which the basis is defined. See Shards documentation http://trilinos.sandia.gov/packages/shards for definition of base cell topology. More... | |
| virtual EBasis | getBasisType () const |
| Returns the basis type. More... | |
| virtual ECoordinates | getCoordinateSystem () const |
| Returns the type of coordinate system for which the basis is defined. More... | |
| virtual int | getDofOrdinal (const int subcDim, const int subcOrd, const int subcDofOrd) |
| DoF tag to ordinal lookup. More... | |
|
virtual const std::vector < std::vector< std::vector < int > > > & | getDofOrdinalData () |
| DoF tag to ordinal data structure. | |
| virtual const std::vector< int > & | getDofTag (const int dofOrd) |
| DoF ordinal to DoF tag lookup. More... | |
| virtual const std::vector < std::vector< int > > & | getAllDofTags () |
| Retrieves all DoF tags. More... | |
Private Member Functions | |
| void | initializeTags () |
| Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays. | |
Additional Inherited Members | |
Protected Attributes inherited from Intrepid::Basis< Scalar, ArrayScalar > | |
| int | basisCardinality_ |
| Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. | |
| int | 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 http://trilinos.sandia.gov/packages/shards for definition of base cell topology. | |
| EBasis | basisType_ |
| Type of the basis. | |
| ECoordinates | basisCoordinates_ |
| The coordinate system for which the basis is defined. | |
| bool | basisTagsAreSet_ |
| "true" if tagToOrdinal_ and ordinalToTag_ have been initialized | |
| std::vector< std::vector< int > > | ordinalToTag_ |
| DoF ordinal to tag lookup table. More... | |
| std::vector< std::vector < std::vector< int > > > | tagToOrdinal_ |
| DoF tag to ordinal lookup table. More... | |
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Wedge cell.
Implements Nedelec basis of degree 1 on the reference Wedge cell. The basis has
cardinality 9 and spans an INCOMPLETE tri-linear 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 | 1 | 0 | 0 | 1 | L_0(u) = (u.t)(0.5,0,-1) | |---------|--------------|--------------|--------------|-------------|----------------------------| | 1 | 1 | 1 | 0 | 1 | L_1(u) = (u.t)(0.5,0.5,-1)| |---------|--------------|--------------|--------------|-------------|----------------------------| | 2 | 1 | 2 | 0 | 1 | L_2(u) = (u.t)(0,0.5,-1) | |---------|--------------|--------------|--------------|-------------|----------------------------| | 3 | 1 | 3 | 0 | 1 | L_3(u) = (u.t)(0.5,0,1) | |---------|--------------|--------------|--------------|-------------|----------------------------| | 4 | 1 | 4 | 0 | 1 | L_4(u) = (u.t)(0.5,0.5,1) | |---------|--------------|--------------|--------------|-------------|----------------------------| | 5 | 1 | 5 | 0 | 1 | L_5(u) = (u.t)(0,0.5,1) | |---------|--------------|--------------|--------------|-------------|----------------------------| | 6 | 1 | 6 | 0 | 1 | L_6(u) = (u.t)(0,0,0) | |---------|--------------|--------------|--------------|-------------|----------------------------| | 7 | 1 | 7 | 0 | 1 | L_7(u) = (u.t)(1,0,0) | |---------|--------------|--------------|--------------|-------------|----------------------------| | 8 | 1 | 8 | 0 | 1 | L_8(u) = (u.t)(0,1,0) | |=========|==============|==============|==============|=============|============================| | MAX | maxScDim=1 | maxScOrd=8 | maxDfOrd=0 | - | | |=========|==============|==============|==============|=============|============================|
is an edge tangent. Direction of edge tangents follows the vertex order of the edges in the cell topology and runs from edge vertex 0 to edge vertex 1, whereas their length is set equal to the edge length. For example, edge 8 of all Wedge reference cells has vertex order {2,5}, i.e., its tangent runs from vertex 2 of the reference Wedge to vertex 5 of that cell. On the reference Wedge the coordinates of these vertices are (0,1,-1) and (0,1,1), respectively. Therefore, the tangent to edge 8 is (0,1,1) - (0,1,-1) = (0, 0, 2). Because its length already equals edge length, no further rescaling of the edge tangent is needed.Definition at line 109 of file Intrepid_HCURL_WEDGE_I1_FEM.hpp.
|
virtual |
Evaluation of a FEM basis on a reference Wedge cell.
Returns values of <var>operatorType</var> acting on FEM basis functions for a set of
points in the <strong>reference Wedge</strong> cell. For rank and dimensions of
I/O array arguments see Section \ref basis_md_array_sec.
| outputValues | [out] - rank-3 array with the computed basis values |
| inputPoints | [in] - rank-2 array with dimensions (P,D) containing reference points |
| operatorType | [in] - operator applied to basis functions |
Implements Intrepid::Basis< Scalar, ArrayScalar >.
Definition at line 97 of file Intrepid_HCURL_WEDGE_I1_FEMDef.hpp.
1.8.5