49 #ifndef INTREPID_HGRAD_TET_Cn_FEM_ORTH_HPP 
   50 #define INTREPID_HGRAD_TET_Cn_FEM_ORTH_HPP 
   66 template<
class Scalar, 
class ArrayScalar> 
 
   91   void getValues(ArrayScalar &          outputValues,
 
   92                  const ArrayScalar &    inputPoints,
 
   93                  const EOperator        operatorType) 
const;
 
   98   void getValues(ArrayScalar &          outputValues,
 
   99                  const ArrayScalar &    inputPoints,
 
  100                  const ArrayScalar &    cellVertices,
 
  101                  const EOperator        operatorType = OPERATOR_VALUE) 
const;
 
  118 template<
typename Scalar,
typename ArrayScalar, 
unsigned derivOrder>
 
  131   static void tabulate( ArrayScalar & outputValues ,
 
  133                         const ArrayScalar &inputPoints );
 
  141 template<
typename Scalar,
typename ArrayScalar>
 
  153   static void tabulate( ArrayScalar & outputValues ,
 
  155                         const ArrayScalar &inputPoints );
 
  162   static int idx(
int p, 
int q,
int r)
 
  164     return (p+q+r)*(p+q+r+1)*(p+q+r+2)/6+(q+r)*(q+r+1)/2+r;
 
  184   static void jrc( 
const Scalar &alpha , 
const Scalar &beta , 
 
  186                    Scalar &an , Scalar &bn, Scalar &cn )
 
  189     an = (2.0 * n + 1.0 + alpha + beta) * ( 2.0 * n + 2.0 + alpha + beta ) 
 
  190       / ( 2.0 * ( n + 1 ) * ( n + 1 + alpha + beta ) );
 
  191     bn = (alpha*alpha-beta*beta)*(2.0*n+1.0+alpha+beta) 
 
  192       / ( 2.0*(n+1.0)*(2.0*n+alpha+beta)*(n+1.0+alpha+beta) );
 
  193     cn = (n+alpha)*(n+beta)*(2.0*n+2.0+alpha+beta) 
 
  194       / ( (n+1.0)*(n+1.0+alpha+beta)*(2.0*n+alpha+beta) );
 
  209 template<
typename Scalar,
typename ArrayScalar>
 
  222   static void tabulate( ArrayScalar & outputValues ,
 
  224                         const ArrayScalar &inputPoints );
 
Basis_HGRAD_TET_Cn_FEM_ORTH(int degree)
Constructor. 
static void jrc(const Scalar &alpha, const Scalar &beta, const int &n, Scalar &an, Scalar &bn, Scalar &cn)
function for computing the Jacobi recurrence coefficients so that 
static int idx(int p, int q, int r)
function for indexing from orthogonal expansion indices into linear space p+q+r = the degree of the p...
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const 
Evaluation of a FEM basis on a reference Tetrahedron cell. 
Definition file for FEM orthogonal basis functions of arbitrary degree for H(grad) functions on TET...
Header file for the abstract base class Intrepid::Basis. 
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
void initializeTags()
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays. 
This is an internal class with a static member function for tabulating derivatives of orthogonal expa...
Implementation of the default H(grad)-compatible orthogonal basis of arbitrary degree on tetrahedron...
static void tabulate(ArrayScalar &outputValues, const int deg, const ArrayScalar &inputPoints)
basic tabulate mathod evaluates the derivOrder^th derivatives of the basis functions at inputPoints i...