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: PetscFunctionBegin;
179: *factorial = -1.0;
180: PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %" PetscInt_FMT, n);
181: for (PetscInt i = 1; i < n + 1; ++i) f *= (PetscReal)i;
182: *factorial = f;
183: PetscFunctionReturn(PETSC_SUCCESS);
184: }
186: /*MC
187: PetscDTFactorialInt - Compute n! as an integer
189: Input Parameter:
190: . n - a non-negative integer
192: Output Parameter:
193: . factorial - n!
195: Level: beginner
197: 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.
198: M*/
199: static inline PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
200: {
201: PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};
203: PetscFunctionBegin;
204: *factorial = -1;
205: 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);
206: if (n <= 12) {
207: *factorial = facLookup[n];
208: } else {
209: PetscInt f = facLookup[12];
210: PetscInt i;
212: for (i = 13; i < n + 1; ++i) f *= i;
213: *factorial = f;
214: }
215: PetscFunctionReturn(PETSC_SUCCESS);
216: }
218: /*MC
219: PetscDTBinomial - Approximate the binomial coefficient "n choose k"
221: Input Parameters:
222: + n - a non-negative integer
223: - k - an integer between 0 and n, inclusive
225: Output Parameter:
226: . binomial - approximation of the binomial coefficient n choose k
228: Level: beginner
229: M*/
230: static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
231: {
232: PetscFunctionBeginHot;
233: *binomial = -1.0;
234: 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);
235: if (n <= 3) {
236: PetscInt binomLookup[4][4] = {
237: {1, 0, 0, 0},
238: {1, 1, 0, 0},
239: {1, 2, 1, 0},
240: {1, 3, 3, 1}
241: };
243: *binomial = (PetscReal)binomLookup[n][k];
244: } else {
245: PetscReal binom = 1.0;
247: k = PetscMin(k, n - k);
248: for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
249: *binomial = binom;
250: }
251: PetscFunctionReturn(PETSC_SUCCESS);
252: }
254: /*MC
255: PetscDTBinomialInt - Compute the binomial coefficient "n choose k"
257: Input Parameters:
258: + n - a non-negative integer
259: - k - an integer between 0 and n, inclusive
261: Output Parameter:
262: . binomial - the binomial coefficient n choose k
264: Note: this is limited by integers that can be represented by `PetscInt`
266: Level: beginner
267: M*/
268: static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
269: {
270: PetscInt bin;
272: PetscFunctionBegin;
273: *binomial = -1;
274: 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);
275: 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);
276: if (n <= 3) {
277: PetscInt binomLookup[4][4] = {
278: {1, 0, 0, 0},
279: {1, 1, 0, 0},
280: {1, 2, 1, 0},
281: {1, 3, 3, 1}
282: };
284: bin = binomLookup[n][k];
285: } else {
286: PetscInt binom = 1;
288: k = PetscMin(k, n - k);
289: for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
290: bin = binom;
291: }
292: *binomial = bin;
293: PetscFunctionReturn(PETSC_SUCCESS);
294: }
296: /*MC
297: PetscDTEnumPerm - Get a permutation of n integers from its encoding into the integers [0, n!) as a sequence of swaps.
299: A permutation can be described by the operations that convert the lists [0, 1, ..., n-1] into the permutation,
300: by a sequence of swaps, where the ith step swaps whatever number is in ith position with a number that is in
301: some position j >= i. This swap is encoded as the difference (j - i). The difference d_i at step i is less than
302: (n - i). This sequence of n-1 differences [d_0, ..., d_{n-2}] is encoded as the number
303: (n-1)! * d_0 + (n-2)! * d_1 + ... + 1! * d_{n-2}.
305: Input Parameters:
306: + n - a non-negative integer (see note about limits below)
307: - k - an integer in [0, n!)
309: Output Parameters:
310: + perm - the permuted list of the integers [0, ..., n-1]
311: - isOdd - if not NULL, returns whether the permutation used an even or odd number of swaps.
313: Note:
314: 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.
316: Level: beginner
317: M*/
318: static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
319: {
320: PetscInt odd = 0;
321: PetscInt i;
322: PetscInt work[PETSC_FACTORIAL_MAX];
323: PetscInt *w;
325: PetscFunctionBegin;
326: if (isOdd) *isOdd = PETSC_FALSE;
327: 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);
328: w = &work[n - 2];
329: for (i = 2; i <= n; i++) {
330: *(w--) = k % i;
331: k /= i;
332: }
333: for (i = 0; i < n; i++) perm[i] = i;
334: for (i = 0; i < n - 1; i++) {
335: PetscInt s = work[i];
336: PetscInt swap = perm[i];
338: perm[i] = perm[i + s];
339: perm[i + s] = swap;
340: odd ^= (!!s);
341: }
342: if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
343: PetscFunctionReturn(PETSC_SUCCESS);
344: }
346: /*MC
347: PetscDTPermIndex - Encode a permutation of n into an integer in [0, n!). This inverts `PetscDTEnumPerm`.
349: Input Parameters:
350: + n - a non-negative integer (see note about limits below)
351: - perm - the permuted list of the integers [0, ..., n-1]
353: Output Parameters:
354: + k - an integer in [0, n!)
355: - isOdd - if not NULL, returns whether the permutation used an even or odd number of swaps.
357: Note:
358: 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.
360: Level: beginner
361: M*/
362: static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
363: {
364: PetscInt odd = 0;
365: PetscInt i, idx;
366: PetscInt work[PETSC_FACTORIAL_MAX];
367: PetscInt iwork[PETSC_FACTORIAL_MAX];
369: PetscFunctionBeginHot;
370: *k = -1;
371: if (isOdd) *isOdd = PETSC_FALSE;
372: 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);
373: for (i = 0; i < n; i++) work[i] = i; /* partial permutation */
374: for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
375: for (idx = 0, i = 0; i < n - 1; i++) {
376: PetscInt j = perm[i];
377: PetscInt icur = work[i];
378: PetscInt jloc = iwork[j];
379: PetscInt diff = jloc - i;
381: idx = idx * (n - i) + diff;
382: /* swap (i, jloc) */
383: work[i] = j;
384: work[jloc] = icur;
385: iwork[j] = i;
386: iwork[icur] = jloc;
387: odd ^= (!!diff);
388: }
389: *k = idx;
390: if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
391: PetscFunctionReturn(PETSC_SUCCESS);
392: }
394: /*MC
395: PetscDTEnumSubset - Get an ordered subset of the integers [0, ..., n - 1] from its encoding as an integers in [0, n choose k).
396: The encoding is in lexicographic order.
398: Input Parameters:
399: + n - a non-negative integer (see note about limits below)
400: . k - an integer in [0, n]
401: - j - an index in [0, n choose k)
403: Output Parameter:
404: . subset - the jth subset of size k of the integers [0, ..., n - 1]
406: Note:
407: Limited by arguments such that n choose k can be represented by `PetscInt`
409: Level: beginner
411: .seealso: `PetscDTSubsetIndex()`
412: M*/
413: static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
414: {
415: PetscInt Nk;
417: PetscFunctionBeginHot;
418: PetscCall(PetscDTBinomialInt(n, k, &Nk));
419: for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
420: PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
421: PetscInt Nminusk = Nk - Nminuskminus;
423: if (j < Nminuskminus) {
424: subset[l++] = i;
425: Nk = Nminuskminus;
426: } else {
427: j -= Nminuskminus;
428: Nk = Nminusk;
429: }
430: }
431: PetscFunctionReturn(PETSC_SUCCESS);
432: }
434: /*MC
435: 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.
436: This is the inverse of `PetscDTEnumSubset`.
438: Input Parameters:
439: + n - a non-negative integer (see note about limits below)
440: . k - an integer in [0, n]
441: - subset - an ordered subset of the integers [0, ..., n - 1]
443: Output Parameter:
444: . index - the rank of the subset in lexicographic order
446: Note:
447: Limited by arguments such that n choose k can be represented by `PetscInt`
449: Level: beginner
451: .seealso: `PetscDTEnumSubset()`
452: M*/
453: static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
454: {
455: PetscInt j = 0, Nk;
457: PetscFunctionBegin;
458: *index = -1;
459: PetscCall(PetscDTBinomialInt(n, k, &Nk));
460: for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
461: PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
462: PetscInt Nminusk = Nk - Nminuskminus;
464: if (subset[l] == i) {
465: l++;
466: Nk = Nminuskminus;
467: } else {
468: j += Nminuskminus;
469: Nk = Nminusk;
470: }
471: }
472: *index = j;
473: PetscFunctionReturn(PETSC_SUCCESS);
474: }
476: /*MC
477: 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.
479: Input Parameters:
480: + n - a non-negative integer (see note about limits below)
481: . k - an integer in [0, n]
482: - j - an index in [0, n choose k)
484: Output Parameters:
485: + perm - the jth subset of size k of the integers [0, ..., n - 1], followed by its complementary set.
486: - isOdd - if not NULL, return whether perm is an even or odd permutation.
488: Note:
489: Limited by arguments such that n choose k can be represented by `PetscInt`
491: Level: beginner
493: .seealso: `PetscDTEnumSubset()`, `PetscDTSubsetIndex()`
494: M*/
495: static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
496: {
497: PetscInt i, l, m, Nk, odd = 0;
498: PetscInt *subcomp = perm + k;
500: PetscFunctionBegin;
501: if (isOdd) *isOdd = PETSC_FALSE;
502: PetscCall(PetscDTBinomialInt(n, k, &Nk));
503: for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
504: PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
505: PetscInt Nminusk = Nk - Nminuskminus;
507: if (j < Nminuskminus) {
508: perm[l++] = i;
509: Nk = Nminuskminus;
510: } else {
511: subcomp[m++] = i;
512: j -= Nminuskminus;
513: odd ^= ((k - l) & 1);
514: Nk = Nminusk;
515: }
516: }
517: for (; i < n; i++) subcomp[m++] = i;
518: if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
519: PetscFunctionReturn(PETSC_SUCCESS);
520: }
522: struct _p_PetscTabulation {
523: PetscInt K; /* Indicates a k-jet, namely tabulated derivatives up to order k */
524: PetscInt Nr; /* The number of tabulation replicas (often 1) */
525: PetscInt Np; /* The number of tabulation points in a replica */
526: PetscInt Nb; /* The number of functions tabulated */
527: PetscInt Nc; /* The number of function components */
528: PetscInt cdim; /* The coordinate dimension */
529: PetscReal **T; /* The tabulation T[K] of functions and their derivatives
530: T[0] = B[Nr*Np][Nb][Nc]: The basis function values at quadrature points
531: T[1] = D[Nr*Np][Nb][Nc][cdim]: The basis function derivatives at quadrature points
532: T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
533: };
534: typedef struct _p_PetscTabulation *PetscTabulation;
536: typedef PetscErrorCode (*PetscProbFunc)(const PetscReal[], const PetscReal[], PetscReal[]);
538: typedef enum {
539: DTPROB_DENSITY_CONSTANT,
540: DTPROB_DENSITY_GAUSSIAN,
541: DTPROB_DENSITY_MAXWELL_BOLTZMANN,
542: DTPROB_NUM_DENSITY
543: } DTProbDensityType;
544: PETSC_EXTERN const char *const DTProbDensityTypes[];
546: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
547: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
548: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
549: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
550: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
551: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
552: PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
553: PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
554: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
555: PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
556: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
557: PETSC_EXTERN PetscErrorCode PetscPDFGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
558: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal[], const PetscReal[], PetscReal[]);
559: PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
560: PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
561: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
562: PETSC_EXTERN PetscErrorCode PetscPDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
563: PETSC_EXTERN PetscErrorCode PetscCDFConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
564: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant2D(const PetscReal[], const PetscReal[], PetscReal[]);
565: PETSC_EXTERN PetscErrorCode PetscPDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
566: PETSC_EXTERN PetscErrorCode PetscCDFConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
567: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant3D(const PetscReal[], const PetscReal[], PetscReal[]);
568: PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);
570: #include <petscvec.h>
572: PETSC_EXTERN PetscErrorCode PetscProbComputeKSStatistic(Vec, PetscProbFunc, PetscReal *);
574: #endif