49 #ifndef INTREPID_HGRAD_TRI_Cn_FEM_ORTHHPP
50 #define INTREPID_HGRAD_TRI_Cn_FEM_ORTHHPP
65 template<
class Scalar,
class ArrayScalar>
90 void getValues(ArrayScalar & outputValues,
91 const ArrayScalar & inputPoints,
92 const EOperator operatorType)
const;
97 void getValues(ArrayScalar & outputValues,
98 const ArrayScalar & inputPoints,
99 const ArrayScalar & cellVertices,
100 const EOperator operatorType = OPERATOR_VALUE)
const;
117 template<
typename Scalar,
typename ArrayScalar,
unsigned derivOrder>
130 static void tabulate( ArrayScalar & outputValues ,
132 const ArrayScalar &inputPoints );
143 template<
typename Scalar,
typename ArrayScalar>
155 static void tabulate( ArrayScalar & outputValues ,
157 const ArrayScalar &inputPoints );
164 static int idx(
int p,
int q)
166 return (p+q)*(p+q+1)/2+q;
187 static void jrc(
const Scalar &alpha ,
const Scalar &beta ,
189 Scalar &an , Scalar &bn, Scalar &cn )
191 an = (2.0 * n + 1.0 + alpha + beta) * ( 2.0 * n + 2.0 + alpha + beta )
192 / ( 2.0 * ( n + 1 ) * ( n + 1 + alpha + beta ) );
193 bn = (alpha*alpha-beta*beta)*(2.0*n+1.0+alpha+beta)
194 / ( 2.0*(n+1.0)*(2.0*n+alpha+beta)*(n+1.0+alpha+beta) );
195 cn = (n+alpha)*(n+beta)*(2.0*n+2.0+alpha+beta)
196 / ( (n+1.0)*(n+1.0+alpha+beta)*(2.0*n+alpha+beta) );
211 template<
typename Scalar,
typename ArrayScalar>
224 static void tabulate( ArrayScalar & outputValues ,
226 const ArrayScalar &inputPoints );
Implementation of the default H(grad)-compatible orthogonal basis (Dubiner) of arbitrary degree on tr...
Basis_HGRAD_TRI_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
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
static int idx(int p, int q)
function for indexing from orthogonal expansion indices into linear space p+q = the degree of the pol...
void initializeTags()
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.
Header file for the abstract base class Intrepid::Basis.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
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...
Definition file for FEM orthogonal basis functions of arbitrary degree for H(grad) functions on TRI...
This is an internal class with a static member function for tabulating derivatives of orthogonal expa...