| Intrepid
    | 
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,1] reference line cell, using Jacobi polynomials. More...
#include <Intrepid_HGRAD_LINE_Cn_FEM_JACOBI.hpp>
 
  
 | Public Member Functions | |
| Basis_HGRAD_LINE_Cn_FEM_JACOBI (int order, Scalar alpha=0, Scalar beta=0) | |
| Constructor. | |
| void | getValues (ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const | 
| Evaluation of a FEM basis on a reference Line 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. | |
| void | setBasisParameters (int n, Scalar alpha=0, Scalar beta=0) | 
| Sets private data basisDegree_, basisCardinality_, jacobiAlpha_, and jacobiBeta_, to n, n+1, alpha, and beta, respectively. | |
|  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. | |
| Private Attributes | |
| Scalar | jacobiAlpha_ | 
| Scalar | jacobiBeta_ | 
| 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 locally H(grad)-compatible FEM basis of variable order on the [-1,1] reference line cell, using Jacobi polynomials.
      Implements Jacobi basis of variable order \form#162 on
      the reference [-1,1] line cell. Jacobi polynomials depend on three parameters
  ,
,  , and
, and  and are defined via the so-called Gamma function by
 and are defined via the so-called Gamma function by 
![\[ P_n^{(\alpha,\beta)} (z) = \frac{\Gamma (\alpha+n+1)}{n!\Gamma (\alpha+\beta+n+1)} \sum_{m=0}^n {n\choose m} \frac{\Gamma (\alpha + \beta + n + m + 1)}{\Gamma (\alpha + m + 1)} \left(\frac{z-1}{2}\right)^m \]](form_167.png) 
 The basis has cardinality  and spans a COMPLETE linear polynomial space. Basis functions are dual to a unisolvent set of degrees of freedom (DoF) enumerated as follows:
 and spans a COMPLETE linear polynomial space. Basis functions are dual to a unisolvent set of degrees of freedom (DoF) enumerated as follows:
| Basis order | DoF tag table | DoF definition | |||
|---|---|---|---|---|---|
| subc dim | subc ordinal | subc DoF tag | subc num DoFs | ||
| 0 | 1 | 0 | 0 | 1 |   | 
| 1 | 1 | 0 | 0-1 | 2 |   | 
| 2 | 1 | 0 | 0-2 | 3 |   | 
| 3 | 1 | 0 | 0-3 | 4 |   | 
| ... | 1 | 0 | ... | ... | ... | 
| n | 1 | 0 | 0-n | n+1 |   | 
For example, for Legendre polynomials ( \form#173), the first 11 bases are given by
| Basis order | DoF tag table | DoF definition | |||
|---|---|---|---|---|---|
| subc dim | subc ordinal | subc DoF tag | subc num DoFs | ||
| 0 | 1 | 0 | 0 | 1 |   | 
| 1 | 1 | 0 | 0-1 | 2 | and:   | 
| 2 | 1 | 0 | 0-2 | 3 | and:   | 
| 3 | 1 | 0 | 0-3 | 4 | and:   | 
| 4 | 1 | 0 | 0-4 | 5 | and:   | 
| 5 | 1 | 0 | 0-5 | 6 | and:   | 
| 6 | 1 | 0 | 0-6 | 7 | and:   | 
| 7 | 1 | 0 | 0-7 | 8 | and:   | 
| 8 | 1 | 0 | 0-8 | 9 | and:   | 
| 9 | 1 | 0 | 0-9 | 10 | and:   | 
| 10 | 1 | 0 | 0-10 | 11 | and:   | 
Definition at line 131 of file Intrepid_HGRAD_LINE_Cn_FEM_JACOBI.hpp.
| 
 | virtual | 
Evaluation of a FEM basis on a reference Line cell.
    Returns values of <var>operatorType</var> acting on FEM basis functions for a set of
    points in the <strong>reference Line</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 69 of file Intrepid_HGRAD_LINE_Cn_FEM_JACOBIDef.hpp.
References Intrepid::IntrepidPolylib::jacobd(), and Intrepid::IntrepidPolylib::jacobfd().
 1.8.5
 1.8.5