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: /* MANSEC = DM */
 11: /* SUBMANSEC = DT */

 13: PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;

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

 18:   Level: beginner

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

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

 27:   Values:
 28: +  `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` - compute the nodes via linear algebra
 29: -  `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON`         - compute the nodes by solving a nonlinear equation with Newton's method

 31:   Level: intermediate

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

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

 44:   Values:
 45: + `PETSCDTNODES_DEFAULT`     - Nodes chosen by PETSc
 46: . `PETSCDTNODES_GAUSSJACOBI` - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
 47: . `PETSCDTNODES_EQUISPACED`  - Nodes equispaced either including the endpoints or excluding them
 48: - `PETSCDTNODES_TANHSINH`    - Nodes at Tanh-Sinh quadrature points

 50:   Level: intermediate

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

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

 66: PETSC_EXTERN const char *const *const PetscDTNodeTypes;

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

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

 88:   Level: intermediate

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

 98: PETSC_EXTERN const char *const *const PetscDTSimplexQuadratureTypes;

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

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

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

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

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

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

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

156:   Level: advanced

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

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

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

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

185: /*MC
186:    PetscDTFactorial - Approximate n! as a real number

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

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

194:    Level: beginner

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

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

210: /*MC
211:    PetscDTFactorialInt - Compute n! as an integer

213:    Input Parameter:
214: .  n - a non-negative integer

216:    Output Parameter:
217: .  factorial - n!

219:    Level: beginner

221:    Note:
222:    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.

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

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

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

245: /*MC
246:    PetscDTBinomial - Approximate the binomial coefficient `n` choose `k`

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

252:    Output Parameter:
253: .  binomial - approximation of the binomial coefficient `n` choose `k`

255:    Level: beginner

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

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

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

283: /*MC
284:    PetscDTBinomialInt - Compute the binomial coefficient `n` choose `k`

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

290:    Output Parameter:
291: .  binomial - the binomial coefficient `n` choose `k`

293:    Level: beginner

295:    Note:
296:    This is limited by integers that can be represented by `PetscInt`.

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

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

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

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

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

330: /* the following inline routines should be not be inline routines and then Fortran binding can be built automatically */
331: #define PeOp

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

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

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

344:    Level: intermediate

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

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

355: .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTPermIndex()`
356: M*/
357: static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PeOp PetscBool *isOdd)
358: {
359:   PetscInt  odd = 0;
360:   PetscInt  i;
361:   PetscInt  work[PETSC_FACTORIAL_MAX];
362:   PetscInt *w;

364:   PetscFunctionBegin;
365:   if (isOdd) *isOdd = PETSC_FALSE;
366:   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);
367:   if (n >= 2) {
368:     w = &work[n - 2];
369:     for (i = 2; i <= n; i++) {
370:       *(w--) = k % i;
371:       k /= i;
372:     }
373:   }
374:   for (i = 0; i < n; i++) perm[i] = i;
375:   for (i = 0; i < n - 1; i++) {
376:     PetscInt s    = work[i];
377:     PetscInt swap = perm[i];

379:     perm[i]     = perm[i + s];
380:     perm[i + s] = swap;
381:     odd ^= (!!s);
382:   }
383:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
384:   PetscFunctionReturn(PETSC_SUCCESS);
385: }

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

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

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

398:    Level: beginner

400:    Note:
401:    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.

403: .seealso: `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`
404: M*/
405: static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PeOp PetscBool *isOdd)
406: {
407:   PetscInt odd = 0;
408:   PetscInt i, idx;
409:   PetscInt work[PETSC_FACTORIAL_MAX];
410:   PetscInt iwork[PETSC_FACTORIAL_MAX];

412:   PetscFunctionBeginHot;
413:   *k = -1;
414:   if (isOdd) *isOdd = PETSC_FALSE;
415:   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);
416:   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
417:   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
418:   for (idx = 0, i = 0; i < n - 1; i++) {
419:     PetscInt j    = perm[i];
420:     PetscInt icur = work[i];
421:     PetscInt jloc = iwork[j];
422:     PetscInt diff = jloc - i;

424:     idx = idx * (n - i) + diff;
425:     /* swap (i, jloc) */
426:     work[i]     = j;
427:     work[jloc]  = icur;
428:     iwork[j]    = i;
429:     iwork[icur] = jloc;
430:     odd ^= (!!diff);
431:   }
432:   *k = idx;
433:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
434:   PetscFunctionReturn(PETSC_SUCCESS);
435: }

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

441:    Input Parameters:
442: +  n - a non-negative integer (see note about limits below)
443: .  k - an integer in [0, n]
444: -  j - an index in [0, n choose k)

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

449:    Level: beginner

451:    Note:
452:    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`

454: .seealso: `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
455: M*/
456: static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
457: {
458:   PetscInt Nk;

460:   PetscFunctionBeginHot;
461:   PetscCall(PetscDTBinomialInt(n, k, &Nk));
462:   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
463:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
464:     PetscInt Nminusk      = Nk - Nminuskminus;

466:     if (j < Nminuskminus) {
467:       subset[l++] = i;
468:       Nk          = Nminuskminus;
469:     } else {
470:       j -= Nminuskminus;
471:       Nk = Nminusk;
472:     }
473:   }
474:   PetscFunctionReturn(PETSC_SUCCESS);
475: }

477: /*MC
478:    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.
479:    This is the inverse of `PetscDTEnumSubset`.

481:    Input Parameters:
482: +  n      - a non-negative integer (see note about limits below)
483: .  k      - an integer in [0, n]
484: -  subset - an ordered subset of the integers [0, ..., n - 1]

486:    Output Parameter:
487: .  index - the rank of the subset in lexicographic order

489:    Level: beginner

491:    Note:
492:    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`

494: .seealso: `PetscDTEnumSubset()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`, `PetscDTPermIndex()`
495: M*/
496: static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
497: {
498:   PetscInt j = 0, Nk;

500:   PetscFunctionBegin;
501:   *index = -1;
502:   PetscCall(PetscDTBinomialInt(n, k, &Nk));
503:   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
504:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
505:     PetscInt Nminusk      = Nk - Nminuskminus;

507:     if (subset[l] == i) {
508:       l++;
509:       Nk = Nminuskminus;
510:     } else {
511:       j += Nminuskminus;
512:       Nk = Nminusk;
513:     }
514:   }
515:   *index = j;
516:   PetscFunctionReturn(PETSC_SUCCESS);
517: }

519: /*MC
520:    PetscDTEnumSplit - 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.

522:    Input Parameters:
523: +  n - a non-negative integer (see note about limits below)
524: .  k - an integer in [0, n]
525: -  j - an index in [0, n choose k)

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

531:    Level: beginner

533:    Note:
534:    Limited by arguments such that `n` choose `k` can be represented by `PetscInt`

536: .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`, `PetscDTFactorial()`, `PetscDTFactorialInt()`, `PetscDTBinomial()`, `PetscDTBinomialInt()`, `PetscDTEnumPerm()`,
537:           `PetscDTPermIndex()`
538: M*/
539: static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PeOp PetscBool *isOdd)
540: {
541:   PetscInt  i, l, m, Nk, odd = 0;
542:   PetscInt *subcomp = PetscSafePointerPlusOffset(perm, k);

544:   PetscFunctionBegin;
545:   if (isOdd) *isOdd = PETSC_FALSE;
546:   PetscCall(PetscDTBinomialInt(n, k, &Nk));
547:   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
548:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
549:     PetscInt Nminusk      = Nk - Nminuskminus;

551:     if (j < Nminuskminus) {
552:       perm[l++] = i;
553:       Nk        = Nminuskminus;
554:     } else {
555:       subcomp[m++] = i;
556:       j -= Nminuskminus;
557:       odd ^= ((k - l) & 1);
558:       Nk = Nminusk;
559:     }
560:   }
561:   for (; i < n; i++) subcomp[m++] = i;
562:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
563:   PetscFunctionReturn(PETSC_SUCCESS);
564: }

566: struct _n_PetscTabulation {
567:   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
568:   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
569:   PetscInt    Np;   /* The number of tabulation points in a replica */
570:   PetscInt    Nb;   /* The number of functions tabulated */
571:   PetscInt    Nc;   /* The number of function components */
572:   PetscInt    cdim; /* The coordinate dimension */
573:   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
574:                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
575:                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
576:                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
577: };

579: /*S
580:    PetscTabulation - PETSc object that manages tabulations for finite element methods.

582:    Level: intermediate

584:    Note:
585:    This is a pointer to a C struct, hence the data in it may be accessed directly.

587:    Fortran Note:
588:    Use `PetscTabulationGetData()` and `PetscTabulationRestoreData()` to access the arrays in the tabulation.

590:    Developer Note:
591:    TODO: put the meaning of the struct fields in this manual page

593: .seealso: `PetscTabulationDestroy()`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()`
594: S*/
595: typedef struct _n_PetscTabulation *PetscTabulation;

597: /*S
598:   PetscProbFn - A prototype of a PDF or CDF used with PETSc probability operations whose names begin with `PetscProb` such as
599:   `PetscProbComputeKSStatistic()`.

601:   Calling Sequence:
602: + x      - input value
603: . scale  - scale factor, I don't know what this is for
604: - result - the value of the PDF or CDF at the input value

606:   Level: beginner

608:   Developer Note:
609:   Why does this take an array argument for `result` when it seems to be able to output a single value?

611: .seealso: `PetscProbComputeKSStatistic()`, `PetscProbComputeKSStatisticWeighted()`, `PetscPDFMaxwellBoltzmann1D()`
612: S*/
613: typedef PetscErrorCode PetscProbFn(const PetscReal x[], const PetscReal scale[], PetscReal result[]);

615: PETSC_EXTERN_TYPEDEF typedef PetscProbFn *PetscProbFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscProbFn*", );

617: typedef enum {
618:   DTPROB_DENSITY_CONSTANT,
619:   DTPROB_DENSITY_GAUSSIAN,
620:   DTPROB_DENSITY_MAXWELL_BOLTZMANN,
621:   DTPROB_NUM_DENSITY
622: } DTProbDensityType;
623: PETSC_EXTERN const char *const DTProbDensityTypes[];

625: PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann1D;
626: PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann1D;
627: PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann2D;
628: PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann2D;
629: PETSC_EXTERN PetscProbFn    PetscPDFMaxwellBoltzmann3D;
630: PETSC_EXTERN PetscProbFn    PetscCDFMaxwellBoltzmann3D;
631: PETSC_EXTERN PetscProbFn    PetscPDFGaussian1D;
632: PETSC_EXTERN PetscProbFn    PetscCDFGaussian1D;
633: PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian1D;
634: PETSC_EXTERN PetscProbFn    PetscPDFGaussian2D;
635: PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian2D;
636: PETSC_EXTERN PetscProbFn    PetscPDFGaussian3D;
637: PETSC_EXTERN PetscProbFn    PetscPDFSampleGaussian3D;
638: PETSC_EXTERN PetscProbFn    PetscPDFConstant1D;
639: PETSC_EXTERN PetscProbFn    PetscCDFConstant1D;
640: PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant1D;
641: PETSC_EXTERN PetscProbFn    PetscPDFConstant2D;
642: PETSC_EXTERN PetscProbFn    PetscCDFConstant2D;
643: PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant2D;
644: PETSC_EXTERN PetscProbFn    PetscPDFConstant3D;
645: PETSC_EXTERN PetscProbFn    PetscCDFConstant3D;
646: PETSC_EXTERN PetscProbFn    PetscPDFSampleConstant3D;
647: PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFn **, PetscProbFn **, PetscProbFn **);

649: #include <petscvec.h>

651: PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFn *, PetscReal *);
652: PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatisticWeighted(Vec, Vec, PetscProbFn *, PetscReal *);
653: PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatisticMagnitude(Vec, PetscProbFn *, PetscReal *);