xref: /petsc/include/petscdt.h (revision dda711d054cd473ec392a02f9c93c457971afd2c)
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 
631a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
641a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
651a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
661a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
671a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
681a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
691a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
70*dda711d0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
711a989b97SToby Isaac PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
721a989b97SToby Isaac 
731a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscInt *binomial)
741a989b97SToby Isaac {
751a989b97SToby Isaac   PetscFunctionBeginHot;
761a989b97SToby 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);
771a989b97SToby Isaac   if (n <= 3) {
781a989b97SToby Isaac     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
791a989b97SToby Isaac 
801a989b97SToby Isaac     *binomial = binomLookup[n][k];
811a989b97SToby Isaac   } else {
821a989b97SToby Isaac     PetscReal binom = 1;
831a989b97SToby Isaac     PetscInt  i;
841a989b97SToby Isaac 
851a989b97SToby Isaac     k = PetscMin(k, n - k);
861a989b97SToby Isaac     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
871a989b97SToby Isaac     *binomial = binom;
881a989b97SToby Isaac   }
891a989b97SToby Isaac   PetscFunctionReturn(0);
901a989b97SToby Isaac }
911a989b97SToby Isaac 
921a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *work, PetscInt *perm, PetscBool *isOdd)
931a989b97SToby Isaac {
941a989b97SToby Isaac   PetscInt  odd = 0;
951a989b97SToby Isaac   PetscInt  i;
961a989b97SToby Isaac   PetscInt *w = &work[n - 2];
971a989b97SToby Isaac 
981a989b97SToby Isaac   PetscFunctionBeginHot;
991a989b97SToby Isaac   for (i = 2; i <= n; i++) {
1001a989b97SToby Isaac     *(w--) = k % i;
1011a989b97SToby Isaac     k /= i;
1021a989b97SToby Isaac   }
1031a989b97SToby Isaac   for (i = 0; i < n; i++) perm[i] = i;
1041a989b97SToby Isaac   for (i = 0; i < n - 1; i++) {
1051a989b97SToby Isaac     PetscInt s = work[i];
1061a989b97SToby Isaac     PetscInt swap = perm[i];
1071a989b97SToby Isaac 
1081a989b97SToby Isaac     perm[i] = perm[i + s];
1091a989b97SToby Isaac     perm[i + s] = swap;
1101a989b97SToby Isaac     odd ^= (!!s);
1111a989b97SToby Isaac   }
1121a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
1131a989b97SToby Isaac   PetscFunctionReturn(0);
1141a989b97SToby Isaac }
1151a989b97SToby Isaac 
1161a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
1171a989b97SToby Isaac {
1181a989b97SToby Isaac   PetscInt       Nk, i, l;
1191a989b97SToby Isaac   PetscErrorCode ierr;
1201a989b97SToby Isaac 
1211a989b97SToby Isaac   PetscFunctionBeginHot;
1221a989b97SToby Isaac   ierr = PetscDTBinomial(n, k, &Nk);CHKERRQ(ierr);
1231a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
1241a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
1251a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
1261a989b97SToby Isaac 
1271a989b97SToby Isaac     if (j < Nminuskminus) {
1281a989b97SToby Isaac       subset[l++] = i;
1291a989b97SToby Isaac       Nk = Nminuskminus;
1301a989b97SToby Isaac     } else {
1311a989b97SToby Isaac       j -= Nminuskminus;
1321a989b97SToby Isaac       Nk = Nminusk;
1331a989b97SToby Isaac     }
1341a989b97SToby Isaac   }
1351a989b97SToby Isaac   PetscFunctionReturn(0);
1361a989b97SToby Isaac }
1371a989b97SToby Isaac 
1381a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
1391a989b97SToby Isaac {
1401a989b97SToby Isaac   PetscInt       i, j = 0, l, Nk;
1411a989b97SToby Isaac   PetscErrorCode ierr;
1421a989b97SToby Isaac 
1431a989b97SToby Isaac   PetscFunctionBeginHot;
1441a989b97SToby Isaac   ierr = PetscDTBinomial(n, k, &Nk);CHKERRQ(ierr);
1451a989b97SToby Isaac   for (i = 0, l = 0; i < n && l < k; i++) {
1461a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
1471a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
1481a989b97SToby Isaac 
1491a989b97SToby Isaac     if (subset[l] == i) {
1501a989b97SToby Isaac       l++;
1511a989b97SToby Isaac       Nk = Nminuskminus;
1521a989b97SToby Isaac     } else {
1531a989b97SToby Isaac       j += Nminuskminus;
1541a989b97SToby Isaac       Nk = Nminusk;
1551a989b97SToby Isaac     }
1561a989b97SToby Isaac   }
1571a989b97SToby Isaac   *index = j;
1581a989b97SToby Isaac   PetscFunctionReturn(0);
1591a989b97SToby Isaac }
1601a989b97SToby Isaac 
1611a989b97SToby Isaac 
1621a989b97SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset, PetscBool *isOdd)
1631a989b97SToby Isaac {
1641a989b97SToby Isaac   PetscInt       i, l, m, *subcomp, Nk;
1651a989b97SToby Isaac   PetscInt       odd;
1661a989b97SToby Isaac   PetscErrorCode ierr;
1671a989b97SToby Isaac 
1681a989b97SToby Isaac   PetscFunctionBeginHot;
1691a989b97SToby Isaac   ierr = PetscDTBinomial(n, k, &Nk);CHKERRQ(ierr);
1701a989b97SToby Isaac   odd = 0;
1711a989b97SToby Isaac   subcomp = &subset[k];
1721a989b97SToby Isaac   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
1731a989b97SToby Isaac     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
1741a989b97SToby Isaac     PetscInt Nminusk = Nk - Nminuskminus;
1751a989b97SToby Isaac 
1761a989b97SToby Isaac     if (j < Nminuskminus) {
1771a989b97SToby Isaac       subset[l++] = i;
1781a989b97SToby Isaac       Nk = Nminuskminus;
1791a989b97SToby Isaac     } else {
1801a989b97SToby Isaac       subcomp[m++] = i;
1811a989b97SToby Isaac       j -= Nminuskminus;
1821a989b97SToby Isaac       odd ^= ((k - l) & 1);
1831a989b97SToby Isaac       Nk = Nminusk;
1841a989b97SToby Isaac     }
1851a989b97SToby Isaac   }
1861a989b97SToby Isaac   for (; i < n; i++) {
1871a989b97SToby Isaac     subcomp[m++] = i;
1881a989b97SToby Isaac   }
1891a989b97SToby Isaac   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
1901a989b97SToby Isaac   PetscFunctionReturn(0);
1911a989b97SToby Isaac }
1921a989b97SToby Isaac 
19337045ce4SJed Brown #endif
194