137045ce4SJed Brown /* 237045ce4SJed Brown Common tools for constructing discretizations 337045ce4SJed Brown */ 426bd1501SBarry Smith #if !defined(PETSCDT_H) 526bd1501SBarry Smith #define PETSCDT_H 637045ce4SJed Brown 737045ce4SJed Brown #include <petscsys.h> 837045ce4SJed Brown 921454ff5SMatthew G. Knepley /*S 1021454ff5SMatthew G. Knepley PetscQuadrature - Quadrature rule for integration. 1121454ff5SMatthew G. Knepley 12329bbf4eSMatthew G. Knepley Level: beginner 1321454ff5SMatthew G. Knepley 1421454ff5SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy() 1521454ff5SMatthew G. Knepley S*/ 1621454ff5SMatthew G. Knepley typedef struct _p_PetscQuadrature *PetscQuadrature; 1721454ff5SMatthew G. Knepley 188272889dSSatish Balay /*E 19916e780bShannah_mairs PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights 208272889dSSatish Balay 218272889dSSatish Balay Level: intermediate 228272889dSSatish Balay 23f2e8fe4dShannah_mairs $ PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra 24d410ae54Shannah_mairs $ PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method 258272889dSSatish Balay 268272889dSSatish Balay E*/ 27f2e8fe4dShannah_mairs typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType; 288272889dSSatish Balay 2921454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *); 30c9638911SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *); 31bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*); 32bcede257SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt); 33a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*); 34a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt); 35a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]); 36a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []); 3721454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer); 3821454ff5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *); 39a0845e3aSMatthew G. Knepley 4089710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *); 4189710940SMatthew G. Knepley 4237045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*); 4337045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*); 44916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*); 45194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*); 46a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 47a6b92713SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*); 4837045ce4SJed Brown 49b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *); 50b3c0f97bSTom Klotz PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *); 51d525116cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *); 52b3c0f97bSTom Klotz 53916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *); 54916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 55916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 56916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 57916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***); 58916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 59916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 60916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 61916e780bShannah_mairs PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***); 62916e780bShannah_mairs 63*1a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 64*1a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 65*1a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 66*1a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *); 67*1a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *); 68*1a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *); 69*1a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *); 70*1a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *); 71*1a989b97SToby Isaac 72*1a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscInt *binomial) 73*1a989b97SToby Isaac { 74*1a989b97SToby Isaac PetscFunctionBeginHot; 75*1a989b97SToby Isaac if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomal arguments (%D %D) must be non-negative, k <= n\n", n, k); 76*1a989b97SToby Isaac if (n <= 3) { 77*1a989b97SToby Isaac PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}}; 78*1a989b97SToby Isaac 79*1a989b97SToby Isaac *binomial = binomLookup[n][k]; 80*1a989b97SToby Isaac } else { 81*1a989b97SToby Isaac PetscReal binom = 1; 82*1a989b97SToby Isaac PetscInt i; 83*1a989b97SToby Isaac 84*1a989b97SToby Isaac k = PetscMin(k, n - k); 85*1a989b97SToby Isaac for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1); 86*1a989b97SToby Isaac *binomial = binom; 87*1a989b97SToby Isaac } 88*1a989b97SToby Isaac PetscFunctionReturn(0); 89*1a989b97SToby Isaac } 90*1a989b97SToby Isaac 91*1a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *work, PetscInt *perm, PetscBool *isOdd) 92*1a989b97SToby Isaac { 93*1a989b97SToby Isaac PetscInt odd = 0; 94*1a989b97SToby Isaac PetscInt i; 95*1a989b97SToby Isaac PetscInt *w = &work[n - 2]; 96*1a989b97SToby Isaac 97*1a989b97SToby Isaac PetscFunctionBeginHot; 98*1a989b97SToby Isaac for (i = 2; i <= n; i++) { 99*1a989b97SToby Isaac *(w--) = k % i; 100*1a989b97SToby Isaac k /= i; 101*1a989b97SToby Isaac } 102*1a989b97SToby Isaac for (i = 0; i < n; i++) perm[i] = i; 103*1a989b97SToby Isaac for (i = 0; i < n - 1; i++) { 104*1a989b97SToby Isaac PetscInt s = work[i]; 105*1a989b97SToby Isaac PetscInt swap = perm[i]; 106*1a989b97SToby Isaac 107*1a989b97SToby Isaac perm[i] = perm[i + s]; 108*1a989b97SToby Isaac perm[i + s] = swap; 109*1a989b97SToby Isaac odd ^= (!!s); 110*1a989b97SToby Isaac } 111*1a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 112*1a989b97SToby Isaac PetscFunctionReturn(0); 113*1a989b97SToby Isaac } 114*1a989b97SToby Isaac 115*1a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset) 116*1a989b97SToby Isaac { 117*1a989b97SToby Isaac PetscInt Nk, i, l; 118*1a989b97SToby Isaac PetscErrorCode ierr; 119*1a989b97SToby Isaac 120*1a989b97SToby Isaac PetscFunctionBeginHot; 121*1a989b97SToby Isaac ierr = PetscDTBinomial(n, k, &Nk);CHKERRQ(ierr); 122*1a989b97SToby Isaac for (i = 0, l = 0; i < n && l < k; i++) { 123*1a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 124*1a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 125*1a989b97SToby Isaac 126*1a989b97SToby Isaac if (j < Nminuskminus) { 127*1a989b97SToby Isaac subset[l++] = i; 128*1a989b97SToby Isaac Nk = Nminuskminus; 129*1a989b97SToby Isaac } else { 130*1a989b97SToby Isaac j -= Nminuskminus; 131*1a989b97SToby Isaac Nk = Nminusk; 132*1a989b97SToby Isaac } 133*1a989b97SToby Isaac } 134*1a989b97SToby Isaac PetscFunctionReturn(0); 135*1a989b97SToby Isaac } 136*1a989b97SToby Isaac 137*1a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index) 138*1a989b97SToby Isaac { 139*1a989b97SToby Isaac PetscInt i, j = 0, l, Nk; 140*1a989b97SToby Isaac PetscErrorCode ierr; 141*1a989b97SToby Isaac 142*1a989b97SToby Isaac PetscFunctionBeginHot; 143*1a989b97SToby Isaac ierr = PetscDTBinomial(n, k, &Nk);CHKERRQ(ierr); 144*1a989b97SToby Isaac for (i = 0, l = 0; i < n && l < k; i++) { 145*1a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 146*1a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 147*1a989b97SToby Isaac 148*1a989b97SToby Isaac if (subset[l] == i) { 149*1a989b97SToby Isaac l++; 150*1a989b97SToby Isaac Nk = Nminuskminus; 151*1a989b97SToby Isaac } else { 152*1a989b97SToby Isaac j += Nminuskminus; 153*1a989b97SToby Isaac Nk = Nminusk; 154*1a989b97SToby Isaac } 155*1a989b97SToby Isaac } 156*1a989b97SToby Isaac *index = j; 157*1a989b97SToby Isaac PetscFunctionReturn(0); 158*1a989b97SToby Isaac } 159*1a989b97SToby Isaac 160*1a989b97SToby Isaac 161*1a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset, PetscBool *isOdd) 162*1a989b97SToby Isaac { 163*1a989b97SToby Isaac PetscInt i, l, m, *subcomp, Nk; 164*1a989b97SToby Isaac PetscInt odd; 165*1a989b97SToby Isaac PetscErrorCode ierr; 166*1a989b97SToby Isaac 167*1a989b97SToby Isaac PetscFunctionBeginHot; 168*1a989b97SToby Isaac ierr = PetscDTBinomial(n, k, &Nk);CHKERRQ(ierr); 169*1a989b97SToby Isaac odd = 0; 170*1a989b97SToby Isaac subcomp = &subset[k]; 171*1a989b97SToby Isaac for (i = 0, l = 0, m = 0; i < n && l < k; i++) { 172*1a989b97SToby Isaac PetscInt Nminuskminus = (Nk * (k - l)) / (n - i); 173*1a989b97SToby Isaac PetscInt Nminusk = Nk - Nminuskminus; 174*1a989b97SToby Isaac 175*1a989b97SToby Isaac if (j < Nminuskminus) { 176*1a989b97SToby Isaac subset[l++] = i; 177*1a989b97SToby Isaac Nk = Nminuskminus; 178*1a989b97SToby Isaac } else { 179*1a989b97SToby Isaac subcomp[m++] = i; 180*1a989b97SToby Isaac j -= Nminuskminus; 181*1a989b97SToby Isaac odd ^= ((k - l) & 1); 182*1a989b97SToby Isaac Nk = Nminusk; 183*1a989b97SToby Isaac } 184*1a989b97SToby Isaac } 185*1a989b97SToby Isaac for (; i < n; i++) { 186*1a989b97SToby Isaac subcomp[m++] = i; 187*1a989b97SToby Isaac } 188*1a989b97SToby Isaac if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE; 189*1a989b97SToby Isaac PetscFunctionReturn(0); 190*1a989b97SToby Isaac } 191*1a989b97SToby Isaac 19237045ce4SJed Brown #endif 193