xref: /petsc/include/petscdt.h (revision 16a9b8de2019c191953d52ed8be9fc0021dcda4c)
1 /*
2   Common tools for constructing discretizations
3 */
4 #pragma once
5 
6 #include <petscsys.h>
7 #include <petscdmtypes.h>
8 #include <petscistypes.h>
9 
10 /* SUBMANSEC = DT */
11 
12 PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;
13 
14 /*S
15   PetscQuadrature - Quadrature rule for numerical integration.
16 
17   Level: beginner
18 
19 .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()`
20 S*/
21 typedef struct _p_PetscQuadrature *PetscQuadrature;
22 
23 /*E
24   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
25 
26   Values:
27 +  `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra
28 -  `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON` - compute the nodes by solving a nonlinear equation with Newton's method
29 
30   Level: intermediate
31 
32 .seealso: `PetscQuadrature`
33 E*/
34 typedef enum {
35   PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,
36   PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
37 } PetscGaussLobattoLegendreCreateType;
38 
39 /*E
40   PetscDTNodeType - A description of strategies for generating nodes (both
41   quadrature nodes and nodes for Lagrange polynomials)
42 
43    Values:
44 +  `PETSCDTNODES_DEFAULT` - Nodes chosen by PETSc
45 .  `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
46 .  `PETSCDTNODES_EQUISPACED` - Nodes equispaced either including the endpoints or excluding them
47 -  `PETSCDTNODES_TANHSINH` - Nodes at Tanh-Sinh quadrature points
48 
49   Level: intermediate
50 
51   Note:
52   A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether
53   the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI`
54   with exponents for the weight function.
55 
56 .seealso: `PetscQuadrature`
57 E*/
58 typedef enum {
59   PETSCDTNODES_DEFAULT = -1,
60   PETSCDTNODES_GAUSSJACOBI,
61   PETSCDTNODES_EQUISPACED,
62   PETSCDTNODES_TANHSINH
63 } PetscDTNodeType;
64 
65 PETSC_EXTERN const char *const *const PetscDTNodeTypes;
66 
67 /*E
68   PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices
69 
70   Values:
71 +  `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc
72 .  `PETSCDTSIMPLEXQUAD_CONIC`   - Quadrature rules constructed as
73                                 conically-warped tensor products of 1D
74                                 Gauss-Jacobi quadrature rules.  These are
75                                 explicitly computable in any dimension for any
76                                 degree, and the tensor-product structure can be
77                                 exploited by sum-factorization methods, but
78                                 they are not efficient in terms of nodes per
79                                 polynomial degree.
80 -  `PETSCDTSIMPLEXQUAD_MINSYM`  - Quadrature rules that are fully symmetric
81                                 (symmetries of the simplex preserve the nodes
82                                 and weights) with minimal (or near minimal)
83                                 number of nodes.  In dimensions higher than 1
84                                 these are not simple to compute, so lookup
85                                 tables are used.
86 
87   Level: intermediate
88 
89 .seealso: `PetscQuadrature`, `PetscDTSimplexQuadrature()`
90 E*/
91 typedef enum {
92   PETSCDTSIMPLEXQUAD_DEFAULT = -1,
93   PETSCDTSIMPLEXQUAD_CONIC   = 0,
94   PETSCDTSIMPLEXQUAD_MINSYM
95 } PetscDTSimplexQuadratureType;
96 
97 PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes;
98 
99 PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
100 PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
101 PETSC_EXTERN PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature, DMPolytopeType *);
102 PETSC_EXTERN PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature, DMPolytopeType);
103 PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *);
104 PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
105 PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *);
106 PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
107 PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *);
108 PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]);
109 PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]);
110 PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
111 PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
112 
113 PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
114 PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
115 PETSC_EXTERN PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature, PetscInt *, IS *[]);
116 
117 PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
118 
119 PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
120 PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *);
121 PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
122 PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
123 PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
124 PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *);
125 PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]);
126 PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *);
127 PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
128 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
129 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *);
130 PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
131 PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
132 PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
133 PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *);
134 PETSC_EXTERN PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType, PetscInt, PetscQuadrature *, PetscQuadrature *);
135 
136 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
137 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
138 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
139 
140 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
141 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
142 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
143 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
144 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
145 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
146 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
147 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
148 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
149 
150 PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
151 PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
152 PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
153 PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
154 PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
155 PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
156 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
157 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
158 PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
159 
160 PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *);
161 PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]);
162 PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *);
163 PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]);
164 
165 #if defined(PETSC_USE_64BIT_INDICES)
166   #define PETSC_FACTORIAL_MAX 20
167   #define PETSC_BINOMIAL_MAX  61
168 #else
169   #define PETSC_FACTORIAL_MAX 12
170   #define PETSC_BINOMIAL_MAX  29
171 #endif
172 
173 /*MC
174    PetscDTFactorial - Approximate n! as a real number
175 
176    Input Parameter:
177 .  n - a non-negative integer
178 
179    Output Parameter:
180 .  factorial - n!
181 
182    Level: beginner
183 
184 .seealso: `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
185 M*/
186 static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
187 {
188   PetscReal f = 1.0;
189 
190   PetscFunctionBegin;
191   *factorial = -1.0;
192   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
193   for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i;
194   *factorial = f;
195   PetscFunctionReturn(PETSC_SUCCESS);
196 }
197 
198 /*MC
199    PetscDTFactorialInt - Compute n! as an integer
200 
201    Input Parameter:
202 .  n - a non-negative integer
203 
204    Output Parameter:
205 .  factorial - n!
206 
207    Level: beginner
208 
209    Note:
210    This is limited to n such that n! can be represented by `PetscInt`, which is 12 if `PetscInt` is a signed 32-bit integer and 20 if `PetscInt` is a signed 64-bit integer.
211 
212 .seealso: `PetscDTFactorial()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
213 M*/
214 static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
215 {
216   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
217 
218   PetscFunctionBegin;
219   *factorial = -1;
220   PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of elements %" PetscInt_FMT " is not in supported range [0,%d]", n, PETSC_FACTORIAL_MAX);
221   if (n <= 12) {
222     *factorial = facLookup[n];
223   } else {
224     PetscInt f = facLookup[12];
225     PetscInt i;
226 
227     for (i = 13; i < n + 1; ++i) f *= i;
228     *factorial = f;
229   }
230   PetscFunctionReturn(PETSC_SUCCESS);
231 }
232 
233 /*MC
234    PetscDTBinomial - Approximate the binomial coefficient "n choose k"
235 
236    Input Parameters:
237 +  n - a non-negative integer
238 -  k - an integer between 0 and n, inclusive
239 
240    Output Parameter:
241 .  binomial - approximation of the binomial coefficient n choose k
242 
243    Level: beginner
244 
245 .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
246 M*/
247 static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
248 {
249   PetscFunctionBeginHot;
250   *binomial = -1.0;
251   PetscCheck(n >= 0 && k >= 0 && k <= n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%" PetscInt_FMT " %" PetscInt_FMT ") must be non-negative, k <= n", n, k);
252   if (n <= 3) {
253     PetscInt binomLookup[4][4] = {
254       {1, 0, 0, 0},
255       {1, 1, 0, 0},
256       {1, 2, 1, 0},
257       {1, 3, 3, 1}
258     };
259 
260     *binomial = (PetscReal)binomLookup[n][k];
261   } else {
262     PetscReal binom = 1.0;
263 
264     k = PetscMin(k, n - k);
265     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
266     *binomial = binom;
267   }
268   PetscFunctionReturn(PETSC_SUCCESS);
269 }
270 
271 /*MC
272    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
273 
274    Input Parameters:
275 +  n - a non-negative integer
276 -  k - an integer between 0 and n, inclusive
277 
278    Output Parameter:
279 .  binomial - the binomial coefficient n choose k
280 
281    Level: beginner
282 
283    Note:
284    This is limited by integers that can be represented by `PetscInt`.
285 
286    Use `PetscDTBinomial()` for real number approximations of larger values
287 
288 .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTEnumPerm()`
289 M*/
290 static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
291 {
292   PetscInt bin;
293 
294   PetscFunctionBegin;
295   *binomial = -1;
296   PetscCheck(n >= 0 && k >= 0 && k <= n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%" PetscInt_FMT " %" PetscInt_FMT ") must be non-negative, k <= n", n, k);
297   PetscCheck(n <= PETSC_BINOMIAL_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial elements %" PetscInt_FMT " is larger than max for PetscInt, %d", n, PETSC_BINOMIAL_MAX);
298   if (n <= 3) {
299     PetscInt binomLookup[4][4] = {
300       {1, 0, 0, 0},
301       {1, 1, 0, 0},
302       {1, 2, 1, 0},
303       {1, 3, 3, 1}
304     };
305 
306     bin = binomLookup[n][k];
307   } else {
308     PetscInt binom = 1;
309 
310     k = PetscMin(k, n - k);
311     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
312     bin = binom;
313   }
314   *binomial = bin;
315   PetscFunctionReturn(PETSC_SUCCESS);
316 }
317 
318 /*MC
319    PetscDTEnumPerm - Get a permutation of `n` integers from its encoding into the integers [0, n!) as a sequence of swaps.
320 
321    Input Parameters:
322 +  n - a non-negative integer (see note about limits below)
323 -  k - an integer in [0, n!)
324 
325    Output Parameters:
326 +  perm - the permuted list of the integers [0, ..., n-1]
327 -  isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps.
328 
329    Level: intermediate
330 
331    Notes:
332    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
333    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
334    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
335    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
336    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
337 
338    Limited to `n` such that `n`! can be represented by `PetscInt`, which is 12 if `PetscInt` is a signed 32-bit integer and 20 if `PetscInt` is a signed 64-bit integer.
339 
340 .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()`
341 M*/
342 static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
343 {
344   PetscInt  odd = 0;
345   PetscInt  i;
346   PetscInt  work[PETSC_FACTORIAL_MAX];
347   PetscInt *w;
348 
349   PetscFunctionBegin;
350   if (isOdd) *isOdd = PETSC_FALSE;
351   PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of elements %" PetscInt_FMT " is not in supported range [0,%d]", n, PETSC_FACTORIAL_MAX);
352   w = &work[n - 2];
353   for (i = 2; i <= n; i++) {
354     *(w--) = k % i;
355     k /= i;
356   }
357   for (i = 0; i < n; i++) perm[i] = i;
358   for (i = 0; i < n - 1; i++) {
359     PetscInt s    = work[i];
360     PetscInt swap = perm[i];
361 
362     perm[i]     = perm[i + s];
363     perm[i + s] = swap;
364     odd ^= (!!s);
365   }
366   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
367   PetscFunctionReturn(PETSC_SUCCESS);
368 }
369 
370 /*MC
371    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts `PetscDTEnumPerm()`.
372 
373    Input Parameters:
374 +  n - a non-negative integer (see note about limits below)
375 -  perm - the permuted list of the integers [0, ..., n-1]
376 
377    Output Parameters:
378 +  k - an integer in [0, n!)
379 -  isOdd - if not `NULL`, returns whether the permutation used an even or odd number of swaps.
380 
381    Level: beginner
382 
383    Note:
384    Limited to `n` such that `n`! can be represented by `PetscInt`, which is 12 if `PetscInt` is a signed 32-bit integer and 20 if `PetscInt` is a signed 64-bit integer.
385 
386 .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
387 M*/
388 static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
389 {
390   PetscInt odd = 0;
391   PetscInt i, idx;
392   PetscInt work[PETSC_FACTORIAL_MAX];
393   PetscInt iwork[PETSC_FACTORIAL_MAX];
394 
395   PetscFunctionBeginHot;
396   *k = -1;
397   if (isOdd) *isOdd = PETSC_FALSE;
398   PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of elements %" PetscInt_FMT " is not in supported range [0,%d]", n, PETSC_FACTORIAL_MAX);
399   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
400   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
401   for (idx = 0, i = 0; i < n - 1; i++) {
402     PetscInt j    = perm[i];
403     PetscInt icur = work[i];
404     PetscInt jloc = iwork[j];
405     PetscInt diff = jloc - i;
406 
407     idx = idx * (n - i) + diff;
408     /* swap (i, jloc) */
409     work[i]     = j;
410     work[jloc]  = icur;
411     iwork[j]    = i;
412     iwork[icur] = jloc;
413     odd ^= (!!diff);
414   }
415   *k = idx;
416   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
417   PetscFunctionReturn(PETSC_SUCCESS);
418 }
419 
420 /*MC
421    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
422    The encoding is in lexicographic order.
423 
424    Input Parameters:
425 +  n - a non-negative integer (see note about limits below)
426 .  k - an integer in [0, n]
427 -  j - an index in [0, n choose k)
428 
429    Output Parameter:
430 .  subset - the jth subset of size k of the integers [0, ..., n - 1]
431 
432    Level: beginner
433 
434    Note:
435    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
436 
437 .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
438 M*/
439 static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
440 {
441   PetscInt Nk;
442 
443   PetscFunctionBeginHot;
444   PetscCall(PetscDTBinomialInt(n, k, &Nk));
445   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
446     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
447     PetscInt Nminusk      = Nk - Nminuskminus;
448 
449     if (j < Nminuskminus) {
450       subset[l++] = i;
451       Nk          = Nminuskminus;
452     } else {
453       j -= Nminuskminus;
454       Nk = Nminusk;
455     }
456   }
457   PetscFunctionReturn(PETSC_SUCCESS);
458 }
459 
460 /*MC
461    PetscDTSubsetIndex - Convert an ordered subset of k integers from the set [0, ..., n - 1] to its encoding as an integers in [0, n choose k) in lexicographic order.
462    This is the inverse of `PetscDTEnumSubset`.
463 
464    Input Parameters:
465 +  n - a non-negative integer (see note about limits below)
466 .  k - an integer in [0, n]
467 -  subset - an ordered subset of the integers [0, ..., n - 1]
468 
469    Output Parameter:
470 .  index - the rank of the subset in lexicographic order
471 
472    Level: beginner
473 
474    Note:
475    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
476 
477 .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
478 M*/
479 static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
480 {
481   PetscInt j = 0, Nk;
482 
483   PetscFunctionBegin;
484   *index = -1;
485   PetscCall(PetscDTBinomialInt(n, k, &Nk));
486   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
487     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
488     PetscInt Nminusk      = Nk - Nminuskminus;
489 
490     if (subset[l] == i) {
491       l++;
492       Nk = Nminuskminus;
493     } else {
494       j += Nminuskminus;
495       Nk = Nminusk;
496     }
497   }
498   *index = j;
499   PetscFunctionReturn(PETSC_SUCCESS);
500 }
501 
502 /*MC
503    PetscDTEnumSubset - Split the integers [0, ..., n - 1] into two complementary ordered subsets, the first subset of size k and being the jth subset of that size in lexicographic order.
504 
505    Input Parameters:
506 +  n - a non-negative integer (see note about limits below)
507 .  k - an integer in [0, n]
508 -  j - an index in [0, n choose k)
509 
510    Output Parameters:
511 +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
512 -  isOdd - if not `NULL`, return whether perm is an even or odd permutation.
513 
514    Level: beginner
515 
516    Note:
517    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`
518 
519 .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`,
520           `PetscDTPermIndex()`
521 M*/
522 static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
523 {
524   PetscInt  i, l, m, Nk, odd = 0;
525   PetscInt *subcomp = perm + k;
526 
527   PetscFunctionBegin;
528   if (isOdd) *isOdd = PETSC_FALSE;
529   PetscCall(PetscDTBinomialInt(n, k, &Nk));
530   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
531     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
532     PetscInt Nminusk      = Nk - Nminuskminus;
533 
534     if (j < Nminuskminus) {
535       perm[l++] = i;
536       Nk        = Nminuskminus;
537     } else {
538       subcomp[m++] = i;
539       j -= Nminuskminus;
540       odd ^= ((k - l) & 1);
541       Nk = Nminusk;
542     }
543   }
544   for (; i < n; i++) subcomp[m++] = i;
545   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
546   PetscFunctionReturn(PETSC_SUCCESS);
547 }
548 
549 struct _p_PetscTabulation {
550   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
551   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
552   PetscInt    Np;   /* The number of tabulation points in a replica */
553   PetscInt    Nb;   /* The number of functions tabulated */
554   PetscInt    Nc;   /* The number of function components */
555   PetscInt    cdim; /* The coordinate dimension */
556   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
557                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
558                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
559                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
560 };
561 typedef struct _p_PetscTabulation *PetscTabulation;
562 
563 typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]);
564 
565 typedef enum {
566   DTPROB_DENSITY_CONSTANT,
567   DTPROB_DENSITY_GAUSSIAN,
568   DTPROB_DENSITY_MAXWELL_BOLTZMANN,
569   DTPROB_NUM_DENSITY
570 } DTProbDensityType;
571 PETSC_EXTERN const char *const DTProbDensityTypes[];
572 
573 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
574 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
575 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
576 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
577 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
578 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
579 PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
580 PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
581 PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
582 PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
583 PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
584 PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
585 PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
586 PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
587 PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
588 PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
589 PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
590 PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
591 PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
592 PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
593 PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
594 PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
595 PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);
596 
597 #include <petscvec.h>
598 
599 PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *);
600