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