49 #ifndef Intrepid2_Polynomials_h 
   50 #define Intrepid2_Polynomials_h 
   88     template<
typename OutputValueViewType, 
typename ScalarType>
 
   89     KOKKOS_INLINE_FUNCTION 
void legendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
 
   91       if (n >= 0) outputValues(0) = 1.0;
 
   92       if (n >= 1) outputValues(1) = x;
 
   93       for (
int i=2; i<=n; i++)
 
   95         const ScalarType i_scalar = ScalarType(i);
 
   96         outputValues(i) = (2. - 1. / i_scalar) * x * outputValues(i-1) - (1. - 1. / i_scalar) * outputValues(i-2);
 
  107     template<
typename OutputValueViewType, 
typename ScalarType>
 
  108     KOKKOS_INLINE_FUNCTION 
void legendreDerivativeValues(OutputValueViewType outputValues, 
const OutputValueViewType legendreValues, Intrepid2::ordinal_type n, ScalarType x)
 
  110       if (n >= 0) outputValues(0) = 0.0;
 
  111       if (n >= 1) outputValues(1) = 1.0;
 
  112       for (
int i=2; i<=n; i++)
 
  114         const ScalarType i_scalar = ScalarType(i);
 
  115         outputValues(i) = outputValues(i-2) + (2. * i_scalar - 1.) * legendreValues(i-1);
 
  130     template<
typename OutputValueViewType, 
typename Po
intScalarType>
 
  131     KOKKOS_INLINE_FUNCTION 
void legendreDerivativeValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, PointScalarType x, Intrepid2::ordinal_type dn)
 
  133       const OutputValueViewType nullOutputScalarView;
 
  134       const double alpha = 0;
 
  135       const double beta = 0;
 
  136       const int numPoints = 1;
 
  138       using Layout = 
typename NaturalLayoutForType<PointScalarType>::layout;
 
  140       using UnmanagedPointScalarView = Kokkos::View<PointScalarType*, Layout, Kokkos::MemoryTraits<Kokkos::Unmanaged> >;
 
  141       UnmanagedPointScalarView pointView = UnmanagedPointScalarView(&x,numPoints);
 
  143       for (
int i=0; i<=n; i++)
 
  145         auto jacobiValue = Kokkos::subview(outputValues,Kokkos::pair<Intrepid2::ordinal_type,Intrepid2::ordinal_type>(i,i+1));
 
  146         jacobiValue(0) = 0.0;
 
  149         double scaleFactor = 1.0;
 
  150         for (
int j=1; j<=dn; j++)
 
  152           scaleFactor *= 0.5 * (j+alpha+beta+i);
 
  155         outputValues(i) = jacobiValue(0) * scaleFactor;
 
  168     template<
typename OutputValueViewType, 
typename ScalarType>
 
  169     KOKKOS_INLINE_FUNCTION 
void shiftedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
 
  171       legendreValues(outputValues, n, 2.*x-1.);
 
  184     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  185     KOKKOS_INLINE_FUNCTION 
void shiftedScaledLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  187       using OutputScalar = 
typename OutputValueViewType::value_type;
 
  188       OutputScalar two_x_minus_t = 2. * x - t;
 
  189       OutputScalar t_squared = t * t;
 
  190       if (n >= 0) outputValues(0) = 1.0;
 
  191       if (n >= 1) outputValues(1) = two_x_minus_t;
 
  192       for (
int i=2; i<=n; i++)
 
  194         const ScalarType one_over_i = 1.0 / ScalarType(i);
 
  195         outputValues(i) = one_over_i * ( (2. *i - 1. ) * two_x_minus_t * outputValues(i-1) - (i - 1.) * t_squared * outputValues(i-2));
 
  212     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  213     KOKKOS_INLINE_FUNCTION 
void shiftedScaledIntegratedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  216       shiftedScaledLegendreValues(outputValues,n,x,t);
 
  218       ScalarType P_i_minus_two, P_i_minus_one;
 
  219       if (n >= 0) P_i_minus_two = outputValues(0);
 
  220       if (n >= 1) P_i_minus_one = outputValues(1);
 
  222       if (n >= 0) outputValues(0) = 1.0;
 
  223       if (n >= 1) outputValues(1) = x;
 
  224       for (
int i=2; i<=n; i++)
 
  226         const ScalarType & P_i = outputValues(i); 
 
  227         const ScalarType i_scalar = ScalarType(i);
 
  228         ScalarType L_i = (P_i - t * t * P_i_minus_two) /( 2. * (2. * i_scalar - 1.));
 
  231         P_i_minus_two = P_i_minus_one;
 
  235         outputValues(i) = L_i;
 
  253     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  254     KOKKOS_INLINE_FUNCTION 
void shiftedScaledIntegratedLegendreValues(OutputValueViewType outputValues, 
const OutputValueViewType shiftedScaledLegendreValues,
 
  255                                                                       Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  258       if (n >= 0) outputValues(0) = 1.0;
 
  259       if (n >= 1) outputValues(1) = x;
 
  260       for (
int i=2; i<=n; i++)
 
  262         const ScalarType & P_i           = shiftedScaledLegendreValues(i); 
 
  263         const ScalarType & P_i_minus_two = shiftedScaledLegendreValues(i-2);
 
  264         const ScalarType i_scalar = ScalarType(i);
 
  265         outputValues(i) = (P_i - t * t * P_i_minus_two) /( 2. * (2. * i_scalar - 1.));
 
  307     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  308     KOKKOS_INLINE_FUNCTION 
void shiftedScaledIntegratedLegendreValues_dx(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  310       if (n >= 0) outputValues(0) = 0.0;
 
  311       if (n >= 1) outputValues(1) = 1.0;
 
  312       if (n >= 2) outputValues(2) = 2. * x - t;
 
  313       for (
int i=2; i<=n-1; i++)
 
  315         const ScalarType one_over_i = 1.0 / ScalarType(i);
 
  316         outputValues(i+1) = one_over_i * (2. * i - 1.) * (2. * x - t) * outputValues(i) -  one_over_i * (i - 1.0) * t * t * outputValues(i-1);
 
  329     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  330     KOKKOS_INLINE_FUNCTION 
void shiftedScaledIntegratedLegendreValues_dt(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  333       shiftedScaledLegendreValues(outputValues, n, x, t);
 
  335       ScalarType P_i_minus_2 = outputValues(0);
 
  336       ScalarType P_i_minus_1 = outputValues(1);
 
  338       if (n >= 0) outputValues(0) = 0.0;
 
  339       if (n >= 1) outputValues(1) = 0.0;
 
  341       for (
int i=2; i<=n; i++)
 
  343         const ScalarType L_i_dt = -0.5 * (P_i_minus_1 + t * P_i_minus_2);
 
  345         P_i_minus_2 = P_i_minus_1;
 
  346         P_i_minus_1 = outputValues(i);
 
  348         outputValues(i) = L_i_dt;
 
  362     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  363     KOKKOS_INLINE_FUNCTION 
void shiftedScaledIntegratedLegendreValues_dt(OutputValueViewType outputValues, 
const OutputValueViewType shiftedScaledLegendreValues,
 
  364                                                                          Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  367       if (n >= 0) outputValues(0) = 0.0;
 
  368       if (n >= 1) outputValues(1) = 0.0;
 
  369       for (
int i=2; i<=n; i++)
 
  371         const ScalarType & P_i_minus_1 = shiftedScaledLegendreValues(i-1); 
 
  372         const ScalarType & P_i_minus_2 = shiftedScaledLegendreValues(i-2);
 
  373         outputValues(i) = -0.5 * (P_i_minus_1 + t * P_i_minus_2);
 
  391     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  392     KOKKOS_INLINE_FUNCTION 
void shiftedScaledJacobiValues(OutputValueViewType outputValues, 
double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  394       ScalarType two_x_minus_t = 2. * x - t;
 
  395       ScalarTypeForScaling alpha_squared_t = alpha * alpha * t;
 
  397       if (n >= 0) outputValues(0) = 1.0;
 
  398       if (n >= 1) outputValues(1) = two_x_minus_t + alpha * x;
 
  400       for (
int i=2; i<=n; i++)
 
  402         const ScalarType & P_i_minus_one = outputValues(i-1); 
 
  403         const ScalarType & P_i_minus_two = outputValues(i-2);
 
  405         double a_i = (2. * i) * (i + alpha) * (2. * i + alpha - 2.);
 
  406         double b_i = 2. * i + alpha - 1.;
 
  407         double c_i = (2. * i + alpha) * (2. * i + alpha - 2.);
 
  408         double d_i = 2. * (i + alpha - 1.) * (i - 1.) * (2. * i + alpha);
 
  410         outputValues(i) = (b_i / a_i) * (c_i * two_x_minus_t + alpha_squared_t) * P_i_minus_one - (d_i / a_i) * t * t * P_i_minus_two;
 
  431     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  432     KOKKOS_INLINE_FUNCTION 
void integratedJacobiValues(OutputValueViewType outputValues, 
const OutputValueViewType jacobiValues,
 
  433                                                        double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  436       if (n >= 0) outputValues(0) = 1.0;
 
  437       if (n >= 1) outputValues(1) = x;
 
  439       ScalarType t_squared = t * t;
 
  440       for (
int i=2; i<=n; i++)
 
  442         const ScalarType & P_i         = jacobiValues(i);   
 
  443         const ScalarType & P_i_minus_1 = jacobiValues(i-1);
 
  444         const ScalarType & P_i_minus_2 = jacobiValues(i-2);
 
  446         double a_i = (i + alpha) / ((2. * i + alpha - 1.) * (2. * i + alpha     ));
 
  447         double b_i = alpha       / ((2. * i + alpha - 2.) * (2. * i + alpha     ));
 
  448         double c_i = (i - 1.)    / ((2. * i + alpha - 2.) * (2. * i + alpha - 1.));
 
  450         outputValues(i) = a_i * P_i + b_i * t * P_i_minus_1 - c_i * t_squared * P_i_minus_2;
 
  471     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  472     KOKKOS_INLINE_FUNCTION 
void integratedJacobiValues(OutputValueViewType outputValues,
 
  473                                                        double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  476       shiftedScaledJacobiValues(outputValues, alpha, n, x, t);
 
  478       ScalarType P_i_minus_2 = outputValues(0);
 
  479       ScalarType P_i_minus_1 = outputValues(1);
 
  481       if (n >= 0) outputValues(0) = 1.0;
 
  482       if (n >= 1) outputValues(1) = x;
 
  484       ScalarType t_squared = t * t;
 
  485       for (
int i=2; i<=n; i++)
 
  487         const ScalarType & P_i = outputValues(i);
 
  489         double a_i = (i + alpha) / ((2. * i + alpha - 1.) * (2. * i + alpha     ));
 
  490         double b_i = alpha       / ((2. * i + alpha - 2.) * (2. * i + alpha     ));
 
  491         double c_i = (i - 1.)    / ((2. * i + alpha - 2.) * (2. * i + alpha - 1.));
 
  493         ScalarType L_i = a_i * P_i + b_i * t * P_i_minus_1 - c_i * t_squared * P_i_minus_2;
 
  495         P_i_minus_2 = P_i_minus_1;
 
  498         outputValues(i) = L_i;
 
  511     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  512     KOKKOS_INLINE_FUNCTION 
void integratedJacobiValues_dx(OutputValueViewType outputValues,
 
  513                                                           double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  517       shiftedScaledJacobiValues(outputValues, alpha, n-1, x, t);
 
  520       ScalarType nextValue = 0.0;
 
  521       ScalarType nextNextValue = 0.0;
 
  522       for (
int i=0; i<=n-1; i++)
 
  524         nextNextValue = outputValues(i);
 
  525         outputValues(i) = nextValue;
 
  526         nextValue = nextNextValue;
 
  528       outputValues(n-1) = nextValue;
 
  541     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  542     KOKKOS_INLINE_FUNCTION 
void integratedJacobiValues_dt(OutputValueViewType outputValues, 
const OutputValueViewType jacobiValues,
 
  543                                                           double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  546       if (n >= 0) outputValues(0) = 0.0;
 
  547       if (n >= 1) outputValues(1) = 0.0;
 
  548       for (
int i=2; i<=n; i++)
 
  550         const ScalarType & P_i_minus_1 = jacobiValues(i-1); 
 
  551         const ScalarType & P_i_minus_2 = jacobiValues(i-2);
 
  552         outputValues(i) = - (i-1.) / (2. * i - 2. + alpha) * (P_i_minus_1 + t * P_i_minus_2);
 
  566     template<
typename OutputValueViewType, 
typename ScalarType, 
typename ScalarTypeForScaling>
 
  567     KOKKOS_INLINE_FUNCTION 
void integratedJacobiValues_dt(OutputValueViewType outputValues,
 
  568                                                           double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
 
  571       shiftedScaledJacobiValues(outputValues, alpha, n, x, t);
 
  573       ScalarType P_i_minus_2 = outputValues(0);
 
  574       ScalarType P_i_minus_1 = outputValues(1);
 
  576       if (n >= 0) outputValues(0) = 0.0;
 
  577       if (n >= 1) outputValues(1) = 0.0;
 
  579       for (
int i=2; i<=n; i++)
 
  581         const ScalarType L_i_dt =  - (i-1.) / (2. * i - 2. + alpha) * (P_i_minus_1 + t * P_i_minus_2);
 
  583         P_i_minus_2 = P_i_minus_1;
 
  584         P_i_minus_1 = outputValues(i);
 
  586         outputValues(i) = L_i_dt;
 
Header file for Intrepid2::Polylib class providing orthogonal polynomial calculus and interpolation...
 
Header function for Intrepid2::Util class and other utility functions. 
 
static KOKKOS_INLINE_FUNCTION void JacobiPolynomial(const ordinal_type np, const zViewType z, polyiViewType poly_in, polydViewType polyd, const ordinal_type n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, . 
 
Contains definitions of custom data types in Intrepid2.