Actual source code: petscdt.h

  1: /*
  2:   Common tools for constructing discretizations
  3: */
  4: #ifndef PETSCDT_H
  5: #define PETSCDT_H

  7: #include <petscsys.h>

  9: /* SUBMANSEC = DT */

 11: PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;

 13: /*S
 14:   PetscQuadrature - Quadrature rule for integration.

 16:   Level: beginner

 18: .seealso: `PetscQuadratureCreate()`, `PetscQuadratureDestroy()`
 19: S*/
 20: typedef struct _p_PetscQuadrature *PetscQuadrature;

 22: /*E
 23:   PetscGaussLobattoLegendreCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights

 25:   Level: intermediate

 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

 30: E*/
 31: typedef enum {
 32:   PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,
 33:   PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
 34: } PetscGaussLobattoLegendreCreateType;

 36: /*E
 37:   PetscDTNodeType - A description of strategies for generating nodes (both
 38:   quadrature nodes and nodes for Lagrange polynomials)

 40:   Level: intermediate

 42: $  `PETSCDTNODES_DEFAULT` - Nodes chosen by PETSc
 43: $  `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
 44: $  `PETSCDTNODES_EQUISPACED` - Nodes equispaced either including the endpoints or excluding them
 45: $  `PETSCDTNODES_TANHSINH` - Nodes at Tanh-Sinh quadrature points

 47:   Note:
 48:   A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether
 49:   the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI`
 50:   with exponents for the weight function.

 52: E*/
 53: typedef enum {
 54:   PETSCDTNODES_DEFAULT = -1,
 55:   PETSCDTNODES_GAUSSJACOBI,
 56:   PETSCDTNODES_EQUISPACED,
 57:   PETSCDTNODES_TANHSINH
 58: } PetscDTNodeType;

 60: PETSC_EXTERN const char *const *const PetscDTNodeTypes;

 62: /*E
 63:   PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices

 65:   Level: intermediate

 67: $  `PETSCDTSIMPLEXQUAD_DEFAULT` - Quadrature rule chosen by PETSc
 68: $  `PETSCDTSIMPLEXQUAD_CONIC`   - Quadrature rules constructed as
 69:                                 conically-warped tensor products of 1D
 70:                                 Gauss-Jacobi quadrature rules.  These are
 71:                                 explicitly computable in any dimension for any
 72:                                 degree, and the tensor-product structure can be
 73:                                 exploited by sum-factorization methods, but
 74:                                 they are not efficient in terms of nodes per
 75:                                 polynomial degree.
 76: $  `PETSCDTSIMPLEXQUAD_MINSYM`  - Quadrature rules that are fully symmetric
 77:                                 (symmetries of the simplex preserve the nodes
 78:                                 and weights) with minimal (or near minimal)
 79:                                 number of nodes.  In dimensions higher than 1
 80:                                 these are not simple to compute, so lookup
 81:                                 tables are used.

 83: .seealso: `PetscDTSimplexQuadrature()`
 84: E*/
 85: typedef enum {
 86:   PETSCDTSIMPLEXQUAD_DEFAULT = -1,
 87:   PETSCDTSIMPLEXQUAD_CONIC   = 0,
 88:   PETSCDTSIMPLEXQUAD_MINSYM
 89: } PetscDTSimplexQuadratureType;

 91: PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes;

 93: PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
 94: PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
 95: PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *);
 96: PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
 97: PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *);
 98: PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
 99: PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *);
100: PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]);
101: PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]);
102: PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
103: PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);

105: PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
106: PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);

108: PETSC_EXTERN PetscErrorCode PetscQuadraturePushForward(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, PetscQuadrature *);

110: PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
111: PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *);
112: PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
113: PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
114: PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
115: PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *);
116: PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]);
117: PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *);
118: PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
119: PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
120: PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *);
121: PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
122: PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
123: PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
124: PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *);

126: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
127: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
128: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);

130: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
131: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
132: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
133: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
134: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
135: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
136: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
137: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
138: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);

140: PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
141: PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
142: PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
143: PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
144: PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
145: PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
146: PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
147: PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
148: PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);

150: PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *);
151: PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]);
152: PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *);
153: PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]);

155: #if defined(PETSC_USE_64BIT_INDICES)
156:   #define PETSC_FACTORIAL_MAX 20
157:   #define PETSC_BINOMIAL_MAX  61
158: #else
159:   #define PETSC_FACTORIAL_MAX 12
160:   #define PETSC_BINOMIAL_MAX  29
161: #endif

163: /*MC
164:    PetscDTFactorial - Approximate n! as a real number

166:    Input Parameter:
167: .  n - a non-negative integer

169:    Output Parameter:
170: .  factorial - n!

172:    Level: beginner
173: M*/
174: static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
175: {
176:   PetscReal f = 1.0;

178:   *factorial = -1.0;
180:   for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i;
181:   *factorial = f;
182:   return 0;
183: }

185: /*MC
186:    PetscDTFactorialInt - Compute n! as an integer

188:    Input Parameter:
189: .  n - a non-negative integer

191:    Output Parameter:
192: .  factorial - n!

194:    Level: beginner

196:    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.
197: M*/
198: static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
199: {
200:   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};

202:   *factorial = -1;
204:   if (n <= 12) {
205:     *factorial = facLookup[n];
206:   } else {
207:     PetscInt f = facLookup[12];
208:     PetscInt i;

210:     for (i = 13; i < n + 1; ++i) f *= i;
211:     *factorial = f;
212:   }
213:   return 0;
214: }

216: /*MC
217:    PetscDTBinomial - Approximate the binomial coefficient "n choose k"

219:    Input Parameters:
220: +  n - a non-negative integer
221: -  k - an integer between 0 and n, inclusive

223:    Output Parameter:
224: .  binomial - approximation of the binomial coefficient n choose k

226:    Level: beginner
227: M*/
228: static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
229: {
231:   *binomial = -1.0;
233:   if (n <= 3) {
234:     PetscInt binomLookup[4][4] = {
235:       {1, 0, 0, 0},
236:       {1, 1, 0, 0},
237:       {1, 2, 1, 0},
238:       {1, 3, 3, 1}
239:     };

241:     *binomial = (PetscReal)binomLookup[n][k];
242:   } else {
243:     PetscReal binom = 1.0;

245:     k = PetscMin(k, n - k);
246:     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
247:     *binomial = binom;
248:   }
249:   return 0;
250: }

252: /*MC
253:    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"

255:    Input Parameters:
256: +  n - a non-negative integer
257: -  k - an integer between 0 and n, inclusive

259:    Output Parameter:
260: .  binomial - the binomial coefficient n choose k

262:    Note: this is limited by integers that can be represented by `PetscInt`

264:    Level: beginner
265: M*/
266: static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
267: {
268:   PetscInt bin;

270:   *binomial = -1;
273:   if (n <= 3) {
274:     PetscInt binomLookup[4][4] = {
275:       {1, 0, 0, 0},
276:       {1, 1, 0, 0},
277:       {1, 2, 1, 0},
278:       {1, 3, 3, 1}
279:     };

281:     bin = binomLookup[n][k];
282:   } else {
283:     PetscInt binom = 1;

285:     k = PetscMin(k, n - k);
286:     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
287:     bin = binom;
288:   }
289:   *binomial = bin;
290:   return 0;
291: }

293: /*MC
294:    PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.

296:    A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
297:    by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
298:    some position j >= i.  This swap is encoded as the difference (j - i).  The difference d_i at step i is less than
299:    (n - i).  This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
300:    (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.

302:    Input Parameters:
303: +  n - a non-negative integer (see note about limits below)
304: -  k - an integer in [0, n!)

306:    Output Parameters:
307: +  perm - the permuted list of the integers [0, ..., n-1]
308: -  isOdd - if not NULL, returns whether the permutation used an even or odd number of swaps.

310:    Note:
311:    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.

313:    Level: beginner
314: M*/
315: static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
316: {
317:   PetscInt  odd = 0;
318:   PetscInt  i;
319:   PetscInt  work[PETSC_FACTORIAL_MAX];
320:   PetscInt *w;

322:   if (isOdd) *isOdd = PETSC_FALSE;
324:   w = &work[n - 2];
325:   for (i = 2; i <= n; i++) {
326:     *(w--) = k % i;
327:     k /= i;
328:   }
329:   for (i = 0; i < n; i++) perm[i] = i;
330:   for (i = 0; i < n - 1; i++) {
331:     PetscInt s    = work[i];
332:     PetscInt swap = perm[i];

334:     perm[i]     = perm[i + s];
335:     perm[i + s] = swap;
336:     odd ^= (!!s);
337:   }
338:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
339:   return 0;
340: }

342: /*MC
343:    PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!).  This inverts `PetscDTEnumPerm`.

345:    Input Parameters:
346: +  n - a non-negative integer (see note about limits below)
347: -  perm - the permuted list of the integers [0, ..., n-1]

349:    Output Parameters:
350: +  k - an integer in [0, n!)
351: -  isOdd - if not NULL, returns whether the permutation used an even or odd number of swaps.

353:    Note:
354:    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.

356:    Level: beginner
357: M*/
358: static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
359: {
360:   PetscInt odd = 0;
361:   PetscInt i, idx;
362:   PetscInt work[PETSC_FACTORIAL_MAX];
363:   PetscInt iwork[PETSC_FACTORIAL_MAX];

366:   *k = -1;
367:   if (isOdd) *isOdd = PETSC_FALSE;
369:   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
370:   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
371:   for (idx = 0, i = 0; i < n - 1; i++) {
372:     PetscInt j    = perm[i];
373:     PetscInt icur = work[i];
374:     PetscInt jloc = iwork[j];
375:     PetscInt diff = jloc - i;

377:     idx = idx * (n - i) + diff;
378:     /* swap (i, jloc) */
379:     work[i]     = j;
380:     work[jloc]  = icur;
381:     iwork[j]    = i;
382:     iwork[icur] = jloc;
383:     odd ^= (!!diff);
384:   }
385:   *k = idx;
386:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
387:   return 0;
388: }

390: /*MC
391:    PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
392:    The encoding is in lexicographic order.

394:    Input Parameters:
395: +  n - a non-negative integer (see note about limits below)
396: .  k - an integer in [0, n]
397: -  j - an index in [0, n choose k)

399:    Output Parameter:
400: .  subset - the jth subset of size k of the integers [0, ..., n - 1]

402:    Note:
403:    Limited by arguments such that n choose k can be represented by `PetscInt`

405:    Level: beginner

407: .seealso: `PetscDTSubsetIndex()`
408: M*/
409: static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
410: {
411:   PetscInt Nk;

414:   PetscDTBinomialInt(n, k, &Nk);
415:   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
416:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
417:     PetscInt Nminusk      = Nk - Nminuskminus;

419:     if (j < Nminuskminus) {
420:       subset[l++] = i;
421:       Nk          = Nminuskminus;
422:     } else {
423:       j -= Nminuskminus;
424:       Nk = Nminusk;
425:     }
426:   }
427:   return 0;
428: }

430: /*MC
431:    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.
432:    This is the inverse of `PetscDTEnumSubset`.

434:    Input Parameters:
435: +  n - a non-negative integer (see note about limits below)
436: .  k - an integer in [0, n]
437: -  subset - an ordered subset of the integers [0, ..., n - 1]

439:    Output Parameter:
440: .  index - the rank of the subset in lexicographic order

442:    Note:
443:    Limited by arguments such that n choose k can be represented by `PetscInt`

445:    Level: beginner

447: .seealso: `PetscDTEnumSubset()`
448: M*/
449: static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
450: {
451:   PetscInt j = 0, Nk;

453:   *index = -1;
454:   PetscDTBinomialInt(n, k, &Nk);
455:   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
456:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
457:     PetscInt Nminusk      = Nk - Nminuskminus;

459:     if (subset[l] == i) {
460:       l++;
461:       Nk = Nminuskminus;
462:     } else {
463:       j += Nminuskminus;
464:       Nk = Nminusk;
465:     }
466:   }
467:   *index = j;
468:   return 0;
469: }

471: /*MC
472:    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.

474:    Input Parameters:
475: +  n - a non-negative integer (see note about limits below)
476: .  k - an integer in [0, n]
477: -  j - an index in [0, n choose k)

479:    Output Parameters:
480: +  perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
481: -  isOdd - if not NULL, return whether perm is an even or odd permutation.

483:    Note:
484:    Limited by arguments such that n choose k can be represented by `PetscInt`

486:    Level: beginner

488: .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`
489: M*/
490: static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
491: {
492:   PetscInt  i, l, m, Nk, odd = 0;
493:   PetscInt *subcomp = perm + k;

495:   if (isOdd) *isOdd = PETSC_FALSE;
496:   PetscDTBinomialInt(n, k, &Nk);
497:   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
498:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
499:     PetscInt Nminusk      = Nk - Nminuskminus;

501:     if (j < Nminuskminus) {
502:       perm[l++] = i;
503:       Nk        = Nminuskminus;
504:     } else {
505:       subcomp[m++] = i;
506:       j -= Nminuskminus;
507:       odd ^= ((k - l) & 1);
508:       Nk = Nminusk;
509:     }
510:   }
511:   for (; i < n; i++) subcomp[m++] = i;
512:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
513:   return 0;
514: }

516: struct _p_PetscTabulation {
517:   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
518:   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
519:   PetscInt    Np;   /* The number of tabulation points in a replica */
520:   PetscInt    Nb;   /* The number of functions tabulated */
521:   PetscInt    Nc;   /* The number of function components */
522:   PetscInt    cdim; /* The coordinate dimension */
523:   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
524:                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
525:                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
526:                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
527: };
528: typedef struct _p_PetscTabulation *PetscTabulation;

530: typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]);

532: typedef enum {
533:   DTPROB_DENSITY_CONSTANT,
534:   DTPROB_DENSITY_GAUSSIAN,
535:   DTPROB_DENSITY_MAXWELL_BOLTZMANN,
536:   DTPROB_NUM_DENSITY
537: } DTProbDensityType;
538: PETSC_EXTERN const char *const DTProbDensityTypes[];

540: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
541: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
542: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
543: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
544: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
545: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
546: PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
547: PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
548: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
549: PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
550: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
551: PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
552: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
553: PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
554: PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
555: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
556: PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
557: PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
558: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
559: PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
560: PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
561: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
562: PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);

564: #include <petscvec.h>

566: PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *);

568: #endif