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