xref: /petsc/include/petscdt.h (revision 91e63d38360eb9bc922f79d792328cc4769c01ac)
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 PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
58 PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
59 PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
60 PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
61 
62 PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
63 PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
64 
65 PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);
66 
67 PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
68 PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal,PetscReal,PetscInt,PetscReal *);
69 PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt,PetscReal,PetscReal,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
70 PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal,PetscReal,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]);
71 PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]);
72 PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt,PetscInt,PetscInt,PetscInt*);
73 PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscInt,PetscReal[]);
74 PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
75 PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
76 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
77 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
78 PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
79 PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
80 PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
81 
82 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
83 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
84 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
85 
86 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
87 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
88 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
89 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
90 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
91 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
92 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
93 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
94 PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
95 
96 PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
97 PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
98 PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
99 PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
100 PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
101 PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
102 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
103 PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
104 PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
105 
106 PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt,PetscInt,const PetscInt[],PetscInt*);
107 PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt,PetscInt,PetscInt,PetscInt[]);
108 PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt,const PetscInt[],PetscInt*);
109 PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt,PetscInt,PetscInt[]);
110 
111 #if defined(PETSC_USE_64BIT_INDICES)
112 #define PETSC_FACTORIAL_MAX 20
113 #define PETSC_BINOMIAL_MAX  61
114 #else
115 #define PETSC_FACTORIAL_MAX 12
116 #define PETSC_BINOMIAL_MAX  29
117 #endif
118 
119 /*MC
120    PetscDTFactorial - Approximate n! as a real number
121 
122    Input Parameter:
123 .  n - a non-negative integer
124 
125    Output Parameter:
126 .  factorial - n!
127 
128    Level: beginner
129 M*/
130 static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
131 {
132   PetscReal f = 1.0;
133   PetscInt  i;
134 
135   PetscFunctionBegin;
136   *factorial = -1.0;
137   PetscCheck(n >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D", n);
138   for (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 %D 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 (%D %D) 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     PetscInt  i;
199 
200     k = PetscMin(k, n - k);
201     for (i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
202     *binomial = binom;
203   }
204   PetscFunctionReturn(0);
205 }
206 
207 /*MC
208    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
209 
210    Input Parameters:
211 +  n - a non-negative integer
212 -  k - an integer between 0 and n, inclusive
213 
214    Output Parameter:
215 .  binomial - the binomial coefficient n choose k
216 
217    Note: this is limited by integers that can be represented by PetscInt
218 
219    Level: beginner
220 M*/
221 static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
222 {
223   PetscInt bin;
224 
225   PetscFunctionBegin;
226   *binomial = -1;
227   PetscCheck(n >= 0 && k >= 0 && k <= n,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n", n, k);
228   PetscCheck(n <= PETSC_BINOMIAL_MAX,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial elements %D is larger than max for PetscInt, %D", n, PETSC_BINOMIAL_MAX);
229   if (n <= 3) {
230     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};
231 
232     bin = binomLookup[n][k];
233   } else {
234     PetscInt  binom = 1;
235     PetscInt  i;
236 
237     k = PetscMin(k, n - k);
238     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
239     bin = binom;
240   }
241   *binomial = bin;
242   PetscFunctionReturn(0);
243 }
244 
245 /*MC
246    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
247 
248    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
249    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
250    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
251    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
252    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
253 
254    Input Parameters:
255 +  n - a non-negative integer (see note about limits below)
256 -  k - an integer in [0, n!)
257 
258    Output Parameters:
259 +  perm - the permuted list of the integers [0, ..., n-1]
260 -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
261 
262    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.
263 
264    Level: beginner
265 M*/
266 static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
267 {
268   PetscInt  odd = 0;
269   PetscInt  i;
270   PetscInt  work[PETSC_FACTORIAL_MAX];
271   PetscInt *w;
272 
273   PetscFunctionBegin;
274   if (isOdd) *isOdd = PETSC_FALSE;
275   PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]",n,PETSC_FACTORIAL_MAX);
276   w = &work[n - 2];
277   for (i = 2; i <= n; i++) {
278     *(w--) = k % i;
279     k /= i;
280   }
281   for (i = 0; i < n; i++) perm[i] = i;
282   for (i = 0; i < n - 1; i++) {
283     PetscInt s = work[i];
284     PetscInt swap = perm[i];
285 
286     perm[i] = perm[i + s];
287     perm[i + s] = swap;
288     odd ^= (!!s);
289   }
290   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
291   PetscFunctionReturn(0);
292 }
293 
294 /*MC
295    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts PetscDTEnumPerm.
296 
297    Input Parameters:
298 +  n - a non-negative integer (see note about limits below)
299 -  perm - the permuted list of the integers [0, ..., n-1]
300 
301    Output Parameters:
302 +  k - an integer in [0, n!)
303 -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.
304 
305    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.
306 
307    Level: beginner
308 M*/
309 static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
310 {
311   PetscInt  odd = 0;
312   PetscInt  i, idx;
313   PetscInt  work[PETSC_FACTORIAL_MAX];
314   PetscInt  iwork[PETSC_FACTORIAL_MAX];
315 
316   PetscFunctionBeginHot;
317   *k = -1;
318   if (isOdd) *isOdd = PETSC_FALSE;
319   PetscCheck(n >= 0 && n <= PETSC_FACTORIAL_MAX,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]",n,PETSC_FACTORIAL_MAX);
320   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
321   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
322   for (idx = 0, i = 0; i < n - 1; i++) {
323     PetscInt j = perm[i];
324     PetscInt icur = work[i];
325     PetscInt jloc = iwork[j];
326     PetscInt diff = jloc - i;
327 
328     idx = idx * (n - i) + diff;
329     /* swap (i, jloc) */
330     work[i] = j;
331     work[jloc] = icur;
332     iwork[j] = i;
333     iwork[icur] = jloc;
334     odd ^= (!!diff);
335   }
336   *k = idx;
337   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
338   PetscFunctionReturn(0);
339 }
340 
341 /*MC
342    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
343    The encoding is in lexicographic order.
344 
345    Input Parameters:
346 +  n - a non-negative integer (see note about limits below)
347 .  k - an integer in [0, n]
348 -  j - an index in [0, n choose k)
349 
350    Output Parameter:
351 .  subset - the jth subset of size k of the integers [0, ..., n - 1]
352 
353    Note: this is limited by arguments such that n choose k can be represented by PetscInt
354 
355    Level: beginner
356 
357 .seealso: PetscDTSubsetIndex()
358 M*/
359 static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
360 {
361   PetscInt       Nk, i, l;
362   PetscErrorCode ierr;
363 
364   PetscFunctionBeginHot;
365   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
366   for (i = 0, l = 0; i < n && l < k; i++) {
367     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
368     PetscInt Nminusk = Nk - Nminuskminus;
369 
370     if (j < Nminuskminus) {
371       subset[l++] = i;
372       Nk = Nminuskminus;
373     } else {
374       j -= Nminuskminus;
375       Nk = Nminusk;
376     }
377   }
378   PetscFunctionReturn(0);
379 }
380 
381 /*MC
382    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.
383 
384    Input Parameters:
385 +  n - a non-negative integer (see note about limits below)
386 .  k - an integer in [0, n]
387 -  subset - an ordered subset of the integers [0, ..., n - 1]
388 
389    Output Parameter:
390 .  index - the rank of the subset in lexicographic order
391 
392    Note: this is limited by arguments such that n choose k can be represented by PetscInt
393 
394    Level: beginner
395 
396 .seealso: PetscDTEnumSubset()
397 M*/
398 static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
399 {
400   PetscInt       i, j = 0, l, Nk;
401   PetscErrorCode ierr;
402 
403   PetscFunctionBegin;
404   *index = -1;
405   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
406   for (i = 0, l = 0; i < n && l < k; i++) {
407     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
408     PetscInt Nminusk = Nk - Nminuskminus;
409 
410     if (subset[l] == i) {
411       l++;
412       Nk = Nminuskminus;
413     } else {
414       j += Nminuskminus;
415       Nk = Nminusk;
416     }
417   }
418   *index = j;
419   PetscFunctionReturn(0);
420 }
421 
422 /*MC
423    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.
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 Parameters:
431 +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
432 -  isOdd - if not NULL, return whether perm is an even or odd permutation.
433 
434    Note: this is limited by arguments such that n choose k can be represented by PetscInt
435 
436    Level: beginner
437 
438 .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
439 M*/
440 static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
441 {
442   PetscInt       i, l, m, *subcomp, Nk;
443   PetscInt       odd;
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   if (isOdd) *isOdd = PETSC_FALSE;
448   ierr = PetscDTBinomialInt(n, k, &Nk);CHKERRQ(ierr);
449   odd = 0;
450   subcomp = &perm[k];
451   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
452     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
453     PetscInt Nminusk = Nk - Nminuskminus;
454 
455     if (j < Nminuskminus) {
456       perm[l++] = i;
457       Nk = Nminuskminus;
458     } else {
459       subcomp[m++] = i;
460       j -= Nminuskminus;
461       odd ^= ((k - l) & 1);
462       Nk = Nminusk;
463     }
464   }
465   for (; i < n; i++) {
466     subcomp[m++] = i;
467   }
468   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
469   PetscFunctionReturn(0);
470 }
471 
472 struct _p_PetscTabulation {
473   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
474   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
475   PetscInt    Np;   /* The number of tabulation points in a replica */
476   PetscInt    Nb;   /* The number of functions tabulated */
477   PetscInt    Nc;   /* The number of function components */
478   PetscInt    cdim; /* The coordinate dimension */
479   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
480                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
481                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
482                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
483 };
484 typedef struct _p_PetscTabulation *PetscTabulation;
485 
486 typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]);
487 
488 typedef enum {DTPROB_DENSITY_CONSTANT, DTPROB_DENSITY_GAUSSIAN, DTPROB_DENSITY_MAXWELL_BOLTZMANN, DTPROB_NUM_DENSITY} DTProbDensityType;
489 PETSC_EXTERN const char * const DTProbDensityTypes[];
490 
491 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
492 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
493 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
494 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
495 PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
496 PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
497 PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
498 PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
499 PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
500 PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
501 PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
502 PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
503 PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
504 PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
505 PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);
506 
507 #include <petscvec.h>
508 
509 PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *);
510 
511 #endif
512