Actual source code: petscdt.h

  1: /*
  2:   Common tools for constructing discretizations
  3: */
  4: #pragma once

  6: #include <petscsys.h>
  7: #include <petscdmtypes.h>
  8: #include <petscistypes.h>

 10: /* SUBMANSEC = DT */

 12: PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;

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

 17:   Level: beginner

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

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

 26:   Values:
 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:   Level: intermediate

 32: .seealso: `PetscQuadrature`
 33: E*/
 34: typedef enum {
 35:   PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,
 36:   PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
 37: } PetscGaussLobattoLegendreCreateType;

 39: /*E
 40:   PetscDTNodeType - A description of strategies for generating nodes (both
 41:   quadrature nodes and nodes for Lagrange polynomials)

 43:   Values:
 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

 49:   Level: intermediate

 51:   Note:
 52:   A `PetscDTNodeType` can be paired with a `PetscBool` to indicate whether
 53:   the nodes include endpoints or not, and in the case of `PETSCDT_GAUSSJACOBI`
 54:   with exponents for the weight function.

 56: .seealso: `PetscQuadrature`
 57: E*/
 58: typedef enum {
 59:   PETSCDTNODES_DEFAULT = -1,
 60:   PETSCDTNODES_GAUSSJACOBI,
 61:   PETSCDTNODES_EQUISPACED,
 62:   PETSCDTNODES_TANHSINH
 63: } PetscDTNodeType;

 65: PETSC_EXTERN const char *const *const PetscDTNodeTypes;

 67: /*E
 68:   PetscDTSimplexQuadratureType - A description of classes of quadrature rules for simplices

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

 87:   Level: intermediate

 89: .seealso: `PetscQuadrature`, `PetscDTSimplexQuadrature()`
 90: E*/
 91: typedef enum {
 92:   PETSCDTSIMPLEXQUAD_DEFAULT = -1,
 93:   PETSCDTSIMPLEXQUAD_CONIC   = 0,
 94:   PETSCDTSIMPLEXQUAD_MINSYM
 95: } PetscDTSimplexQuadratureType;

 97: PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes;

 99: PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
100: PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
101: PETSC_EXTERN PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature, DMPolytopeType *);
102: PETSC_EXTERN PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature, DMPolytopeType);
103: PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt *);
104: PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
105: PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt *);
106: PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
107: PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool *);
108: PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt *, PetscInt *, PetscInt *, const PetscReal *[], const PetscReal *[]);
109: PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[]);
110: PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
111: PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);

113: PETSC_EXTERN PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature, PetscQuadrature, PetscQuadrature *);
114: PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
115: PETSC_EXTERN PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature, PetscInt *, IS *[]);

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

119: PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
120: PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal, PetscReal, PetscInt, PetscReal *);
121: PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt, PetscReal, PetscReal, const PetscReal *, PetscInt, const PetscInt *, PetscReal *, PetscReal *, PetscReal *);
122: PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal, PetscReal, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
123: PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscReal[]);
124: PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt, PetscInt, PetscInt, PetscInt *);
125: PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt, PetscInt, const PetscReal[], PetscInt, PetscInt, PetscInt, PetscReal[]);
126: PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt, PetscReal, PetscReal, PetscReal *, PetscReal *);
127: PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
128: PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal *);
129: PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt, PetscGaussLobattoLegendreCreateType, PetscReal *, PetscReal *);
130: PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
131: PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
132: PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt, PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
133: PETSC_EXTERN PetscErrorCode PetscDTSimplexQuadrature(PetscInt, PetscInt, PetscDTSimplexQuadratureType, PetscQuadrature *);
134: PETSC_EXTERN PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType, PetscInt, PetscQuadrature *, PetscQuadrature *);

136: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
137: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
138: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);

140: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
141: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
142: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
143: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
144: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
145: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
146: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
147: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
148: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);

150: /*MC
151:   PETSC_FORM_DEGREE_UNDEFINED - Indicates that a field does not have
152:   a well-defined form degree in exterior calculus.

154:   Level: advanced

156: .seealso: `PetscDTAltV`, `PetscDualSpaceGetFormDegree()`
157: M*/
158: #define PETSC_FORM_DEGREE_UNDEFINED PETSC_INT_MIN

160: PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
161: PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
162: PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
163: PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
164: PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
165: PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
166: PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
167: PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
168: PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);

170: PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt, PetscInt, const PetscInt[], PetscInt *);
171: PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt, PetscInt, PetscInt, PetscInt[]);
172: PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt, const PetscInt[], PetscInt *);
173: PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt, PetscInt, PetscInt[]);

175: #if defined(PETSC_USE_64BIT_INDICES)
176:   #define PETSC_FACTORIAL_MAX 20
177:   #define PETSC_BINOMIAL_MAX  61
178: #else
179:   #define PETSC_FACTORIAL_MAX 12
180:   #define PETSC_BINOMIAL_MAX  29
181: #endif

183: /*MC
184:    PetscDTFactorial - Approximate n! as a real number

186:    Input Parameter:
187: .  n - a non-negative integer

189:    Output Parameter:
190: .  factorial - n!

192:    Level: beginner

194: .seealso: `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
195: M*/
196: static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
197: {
198:   PetscReal f = 1.0;

200:   PetscFunctionBegin;
201:   *factorial = -1.0;
202:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
203:   for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i;
204:   *factorial = f;
205:   PetscFunctionReturn(PETSC_SUCCESS);
206: }

208: /*MC
209:    PetscDTFactorialInt - Compute n! as an integer

211:    Input Parameter:
212: .  n - a non-negative integer

214:    Output Parameter:
215: .  factorial - n!

217:    Level: beginner

219:    Note:
220:    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.

222: .seealso: `PetscDTFactorial()`, `PetscDTBinomialInt()`, `PetscDTBinomial()`
223: M*/
224: static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
225: {
226:   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};

228:   PetscFunctionBegin;
229:   *factorial = -1;
230:   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);
231:   if (n <= 12) {
232:     *factorial = facLookup[n];
233:   } else {
234:     PetscInt f = facLookup[12];
235:     PetscInt i;

237:     for (i = 13; i < n + 1; ++i) f *= i;
238:     *factorial = f;
239:   }
240:   PetscFunctionReturn(PETSC_SUCCESS);
241: }

243: /*MC
244:    PetscDTBinomial - Approximate the binomial coefficient `n` choose `k`

246:    Input Parameters:
247: +  n - a non-negative integer
248: -  k - an integer between 0 and `n`, inclusive

250:    Output Parameter:
251: .  binomial - approximation of the binomial coefficient `n` choose `k`

253:    Level: beginner

255: .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
256: M*/
257: static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
258: {
259:   PetscFunctionBeginHot;
260:   *binomial = -1.0;
261:   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);
262:   if (n <= 3) {
263:     PetscInt binomLookup[4][4] = {
264:       {1, 0, 0, 0},
265:       {1, 1, 0, 0},
266:       {1, 2, 1, 0},
267:       {1, 3, 3, 1}
268:     };

270:     *binomial = (PetscReal)binomLookup[n][k];
271:   } else {
272:     PetscReal binom = 1.0;

274:     k = PetscMin(k, n - k);
275:     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
276:     *binomial = binom;
277:   }
278:   PetscFunctionReturn(PETSC_SUCCESS);
279: }

281: /*MC
282:    PetscDTBinomialInt - Compute the binomial coefficient `n` choose `k`

284:    Input Parameters:
285: +  n - a non-negative integer
286: -  k - an integer between 0 and `n`, inclusive

288:    Output Parameter:
289: .  binomial - the binomial coefficient `n` choose `k`

291:    Level: beginner

293:    Note:
294:    This is limited by integers that can be represented by `PetscInt`.

296:    Use `PetscDTBinomial()` for real number approximations of larger values

298: .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTEnumPerm()`
299: M*/
300: static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
301: {
302:   PetscInt bin;

304:   PetscFunctionBegin;
305:   *binomial = -1;
306:   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);
307:   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);
308:   if (n <= 3) {
309:     PetscInt binomLookup[4][4] = {
310:       {1, 0, 0, 0},
311:       {1, 1, 0, 0},
312:       {1, 2, 1, 0},
313:       {1, 3, 3, 1}
314:     };

316:     bin = binomLookup[n][k];
317:   } else {
318:     PetscInt binom = 1;

320:     k = PetscMin(k, n - k);
321:     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
322:     bin = binom;
323:   }
324:   *binomial = bin;
325:   PetscFunctionReturn(PETSC_SUCCESS);
326: }

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

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

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

339:    Level: intermediate

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

348:    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.

350: .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()`
351: M*/
352: static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
353: {
354:   PetscInt  odd = 0;
355:   PetscInt  i;
356:   PetscInt  work[PETSC_FACTORIAL_MAX];
357:   PetscInt *w;

359:   PetscFunctionBegin;
360:   if (isOdd) *isOdd = PETSC_FALSE;
361:   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);
362:   w = &work[n - 2];
363:   for (i = 2; i <= n; i++) {
364:     *(w--) = k % i;
365:     k /= i;
366:   }
367:   for (i = 0; i < n; i++) perm[i] = i;
368:   for (i = 0; i < n - 1; i++) {
369:     PetscInt s    = work[i];
370:     PetscInt swap = perm[i];

372:     perm[i]     = perm[i + s];
373:     perm[i + s] = swap;
374:     odd ^= (!!s);
375:   }
376:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
377:   PetscFunctionReturn(PETSC_SUCCESS);
378: }

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

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

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

391:    Level: beginner

393:    Note:
394:    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.

396: .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
397: M*/
398: static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
399: {
400:   PetscInt odd = 0;
401:   PetscInt i, idx;
402:   PetscInt work[PETSC_FACTORIAL_MAX];
403:   PetscInt iwork[PETSC_FACTORIAL_MAX];

405:   PetscFunctionBeginHot;
406:   *k = -1;
407:   if (isOdd) *isOdd = PETSC_FALSE;
408:   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);
409:   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
410:   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
411:   for (idx = 0, i = 0; i < n - 1; i++) {
412:     PetscInt j    = perm[i];
413:     PetscInt icur = work[i];
414:     PetscInt jloc = iwork[j];
415:     PetscInt diff = jloc - i;

417:     idx = idx * (n - i) + diff;
418:     /* swap (i, jloc) */
419:     work[i]     = j;
420:     work[jloc]  = icur;
421:     iwork[j]    = i;
422:     iwork[icur] = jloc;
423:     odd ^= (!!diff);
424:   }
425:   *k = idx;
426:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
427:   PetscFunctionReturn(PETSC_SUCCESS);
428: }

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

434:    Input Parameters:
435: +  n - a non-negative integer (see note about limits below)
436: .  k - an integer in [0, n]
437: -  j - an index in [0, n choose k)

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

442:    Level: beginner

444:    Note:
445:    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`

447: .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
448: M*/
449: static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
450: {
451:   PetscInt Nk;

453:   PetscFunctionBeginHot;
454:   PetscCall(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 (j < Nminuskminus) {
460:       subset[l++] = i;
461:       Nk          = Nminuskminus;
462:     } else {
463:       j -= Nminuskminus;
464:       Nk = Nminusk;
465:     }
466:   }
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: /*MC
471:    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.
472:    This is the inverse of `PetscDTEnumSubset`.

474:    Input Parameters:
475: +  n      - a non-negative integer (see note about limits below)
476: .  k      - an integer in [0, n]
477: -  subset - an ordered subset of the integers [0, ..., n - 1]

479:    Output Parameter:
480: .  index - the rank of the subset in lexicographic order

482:    Level: beginner

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

487: .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
488: M*/
489: static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
490: {
491:   PetscInt j = 0, Nk;

493:   PetscFunctionBegin;
494:   *index = -1;
495:   PetscCall(PetscDTBinomialInt(n, k, &Nk));
496:   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
497:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
498:     PetscInt Nminusk      = Nk - Nminuskminus;

500:     if (subset[l] == i) {
501:       l++;
502:       Nk = Nminuskminus;
503:     } else {
504:       j += Nminuskminus;
505:       Nk = Nminusk;
506:     }
507:   }
508:   *index = j;
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }

512: /*MC
513:    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.

515:    Input Parameters:
516: +  n - a non-negative integer (see note about limits below)
517: .  k - an integer in [0, n]
518: -  j - an index in [0, n choose k)

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

524:    Level: beginner

526:    Note:
527:    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`

529: .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`,
530:           `PetscDTPermIndex()`
531: M*/
532: static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
533: {
534:   PetscInt  i, l, m, Nk, odd = 0;
535:   PetscInt *subcomp = perm + k;

537:   PetscFunctionBegin;
538:   if (isOdd) *isOdd = PETSC_FALSE;
539:   PetscCall(PetscDTBinomialInt(n, k, &Nk));
540:   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
541:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
542:     PetscInt Nminusk      = Nk - Nminuskminus;

544:     if (j < Nminuskminus) {
545:       perm[l++] = i;
546:       Nk        = Nminuskminus;
547:     } else {
548:       subcomp[m++] = i;
549:       j -= Nminuskminus;
550:       odd ^= ((k - l) & 1);
551:       Nk = Nminusk;
552:     }
553:   }
554:   for (; i < n; i++) subcomp[m++] = i;
555:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
556:   PetscFunctionReturn(PETSC_SUCCESS);
557: }

559: struct _p_PetscTabulation {
560:   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
561:   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
562:   PetscInt    Np;   /* The number of tabulation points in a replica */
563:   PetscInt    Nb;   /* The number of functions tabulated */
564:   PetscInt    Nc;   /* The number of function components */
565:   PetscInt    cdim; /* The coordinate dimension */
566:   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
567:                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
568:                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
569:                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
570: };
571: typedef struct _p_PetscTabulation *PetscTabulation;

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

575: typedef enum {
576:   DTPROB_DENSITY_CONSTANT,
577:   DTPROB_DENSITY_GAUSSIAN,
578:   DTPROB_DENSITY_MAXWELL_BOLTZMANN,
579:   DTPROB_NUM_DENSITY
580: } DTProbDensityType;
581: PETSC_EXTERN const char *const DTProbDensityTypes[];

583: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
584: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
585: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
586: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
587: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
588: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
589: PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
590: PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
591: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
592: PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
593: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
594: PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
595: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
596: PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
597: PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
598: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
599: PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
600: PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
601: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
602: PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
603: PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
604: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
605: PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);

607: #include <petscvec.h>

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