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