Intrepid
Public Member Functions | Private Member Functions | Private Attributes | List of all members
Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar > Class Template Reference

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>

Inheritance diagram for Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar >:
Intrepid::Basis< Scalar, ArrayScalar >

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...
 

Detailed Description

template<class Scalar, class ArrayScalar>
class Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar >

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#289 on
      the reference [-1,1] line cell. Jacobi polynomials depend on three parameters

$ \alpha $, $ \beta $, and $ n $ 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 \]

The basis has cardinality $n+1$ 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 $ P_0^{(\alpha,\beta)} $
1 1 0 0-1 2 $ P_0^{(\alpha,\beta)}, P_1^{(\alpha,\beta)} $
2 1 0 0-2 3 $ P_0^{(\alpha,\beta)}, P_1^{(\alpha,\beta)}, P_2^{(\alpha,\beta)} $
3 1 0 0-3 4 $ P_0^{(\alpha,\beta)}, P_1^{(\alpha,\beta)}, ..., P_3^{(\alpha,\beta)} $
... 1 0 ... ... ...
n 1 0 0-n n+1 $ P_0^{(\alpha,\beta)}, P_1^{(\alpha,\beta)}, ..., P_n^{(\alpha,\beta)} $
      For example, for Legendre polynomials ( \form#300), 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 1 0 0-1 2 and: $ x $
2 1 0 0-2 3 and: $ \frac{1}{2} (3x^2-1) $
3 1 0 0-3 4 and: $ \frac{1}{2} (5x^3-3x) $
4 1 0 0-4 5 and: $ \frac{1}{8} (35x^4-30x^2+3) $
5 1 0 0-5 6 and: $ \frac{1}{8} (63x^5-70x^3+15x) $
6 1 0 0-6 7 and: $ \frac{1}{16} (231x^6-315x^4+105x^2-5) $
7 1 0 0-7 8 and: $ \frac{1}{16} (429x^7-693x^5+315x^3-35x) $
8 1 0 0-8 9 and: $ \frac{1}{128} (6435x^8-12012x^6+6930x^4-1260x^2+35) $
9 1 0 0-9 10 and: $ \frac{1}{128} (12155x^9-25740x^7+18018x^5-4620x^3+315x) $
10 1 0 0-1011 and: $ \frac{1}{128} (46189x^{10}-109395x^8+90090x^6-30030x^4+3465x^2-63) $

Definition at line 131 of file Intrepid_HGRAD_LINE_Cn_FEM_JACOBI.hpp.

Member Function Documentation

template<class Scalar , class ArrayScalar>
void Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar >::getValues ( ArrayScalar &  outputValues,
const ArrayScalar &  inputPoints,
const EOperator  operatorType 
) const
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 .
Parameters
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().


The documentation for this class was generated from the following files: