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