|
Intrepid
|
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell. More...
#include <Intrepid_HGRAD_TRI_Cn_FEM.hpp>
Public Member Functions | |
| Basis_HGRAD_TRI_Cn_FEM (const int n, const EPointType pointType) | |
| Constructor. | |
| void | getValues (ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const |
| Evaluation of a FEM basis on a reference Triangle cell. | |
| 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. | |
| virtual int | getDegree () const |
| Returns the degree of the basis. | |
| 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. | |
| virtual EBasis | getBasisType () const |
| Returns the basis type. | |
| virtual ECoordinates | getCoordinateSystem () const |
| Returns the type of coordinate system for which the basis is defined. | |
| virtual int | getDofOrdinal (const int subcDim, const int subcOrd, const int subcDofOrd) |
| DoF tag to ordinal lookup. | |
| 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. | |
| virtual const std::vector< std::vector< int > > & | getAllDofTags () |
| Retrieves all DoF tags. | |
Private Member Functions | |
| virtual void | initializeTags () |
| Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays. | |
Private Attributes | |
| Basis_HGRAD_TRI_Cn_FEM_ORTH< Scalar, FieldContainer< Scalar > > | Phis |
| The orthogonal basis on triangles, out of which the nodal basis is constructed. | |
| FieldContainer< Scalar > | V |
| The Vandermonde matrix with V_{ij} = phi_i(x_j), where x_j is the j_th point in the lattice. | |
| FieldContainer< Scalar > | Vinv |
| The inverse of V. The columns of Vinv express the Lagrange basis in terms of the orthogonal basis. | |
| FieldContainer< Scalar > | latticePts |
| stores the points at which degrees of freedom are located. | |
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. | |
| std::vector< std::vector< std::vector< int > > > | tagToOrdinal_ |
| DoF tag to ordinal lookup table. | |
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell.
Implements Lagrangian basis of degree n on the reference Triangle cell. The basis has
cardinality (n+1)(n+2)/2 and spans a COMPLETE polynomial space of degree n.
Basis functions are dual
to a unisolvent set of degrees-of-freedom (DoF) defined on a lattice of order n (see PointTools).
In particular, the degrees of freedom are point evaluation at
\li the vertices
\li (n-1) points on each edge of the triangle
\li max((n-1)(n-2)/2,0) points on the inside of the triangle.
The distribution of these points is specified by the pointType argument to the class constructor.
Currently, either equispaced lattice points or Warburton's warp-blend points are available.
The dof are enumerated according to the ordering on the lattice (see PointTools). In particular,
dof number 0 is at the bottom left vertex (0,0). The dof increase
along the lattice with points along the lines of constant
x adjacent in the enumeration.
Definition at line 84 of file Intrepid_HGRAD_TRI_Cn_FEM.hpp.
| Intrepid::Basis_HGRAD_TRI_Cn_FEM< Scalar, ArrayScalar >::Basis_HGRAD_TRI_Cn_FEM | ( | const int | n, |
| const EPointType | pointType | ||
| ) |
Constructor.
Definition at line 54 of file Intrepid_HGRAD_TRI_Cn_FEMDef.hpp.
|
virtual |
FVD basis evaluation: invocation of this method throws an exception.
Implements Intrepid::Basis< Scalar, ArrayScalar >.
Definition at line 286 of file Intrepid_HGRAD_TRI_Cn_FEMDef.hpp.
|
virtual |
Evaluation of a FEM basis on a reference Triangle cell.
Returns values of <var>operatorType</var> acting on FEM basis functions for a set of
points in the <strong>reference Triangle</strong> cell. For rank and dimensions of
I/O array arguments see Section \ref basis_md_array_sec .
| 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 |
Implements Intrepid::Basis< Scalar, ArrayScalar >.
Definition at line 190 of file Intrepid_HGRAD_TRI_Cn_FEMDef.hpp.
Referenced by main().
|
privatevirtual |
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.
Implements Intrepid::Basis< Scalar, ArrayScalar >.
Definition at line 111 of file Intrepid_HGRAD_TRI_Cn_FEMDef.hpp.
|
private |
stores the points at which degrees of freedom are located.
Definition at line 103 of file Intrepid_HGRAD_TRI_Cn_FEM.hpp.
|
private |
The orthogonal basis on triangles, out of which the nodal basis is constructed.
Definition at line 93 of file Intrepid_HGRAD_TRI_Cn_FEM.hpp.
|
private |
The Vandermonde matrix with V_{ij} = phi_i(x_j), where x_j is the j_th point in the lattice.
Definition at line 96 of file Intrepid_HGRAD_TRI_Cn_FEM.hpp.
|
private |
The inverse of V. The columns of Vinv express the Lagrange basis in terms of the orthogonal basis.
Definition at line 100 of file Intrepid_HGRAD_TRI_Cn_FEM.hpp.