Actual source code: petscdt.h

  1: /*
  2:   Common tools for constructing discretizations
  3: */
  4: #if !defined(PETSCDT_H)
  5: #define PETSCDT_H

  7: #include <petscsys.h>

  9: PETSC_EXTERN PetscClassId PETSCQUADRATURE_CLASSID;

 11: /*S
 12:   PetscQuadrature - Quadrature rule for integration.

 14:   Level: beginner

 16: .seealso:  PetscQuadratureCreate(), PetscQuadratureDestroy()
 17: S*/
 18: typedef struct _p_PetscQuadrature *PetscQuadrature;

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

 23:   Level: intermediate

 25: $  PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA - compute the nodes via linear algebra
 26: $  PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method

 28: E*/
 29: typedef enum {PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON} PetscGaussLobattoLegendreCreateType;

 31: /*E
 32:   PetscDTNodeType - A description of strategies for generating nodes (both
 33:   quadrature nodes and nodes for Lagrange polynomials)

 35:   Level: intermediate

 37: $  PETSCDTNODES_DEFAULT - Nodes chosen by PETSc
 38: $  PETSCDTNODES_GAUSSJACOBI - Nodes at either Gauss-Jacobi or Gauss-Lobatto-Jacobi quadrature points
 39: $  PETSCDTNODES_EQUISPACED - Nodes equispaced either including the endpoints or excluding them
 40: $  PETSCDTNODES_TANHSINH - Nodes at Tanh-Sinh quadrature points

 42:   Note: a PetscDTNodeType can be paired with a PetscBool to indicate whether
 43:   the nodes include endpoints or not, and in the case of PETSCDT_GAUSSJACOBI
 44:   with exponents for the weight function.

 46: E*/
 47: typedef enum {PETSCDTNODES_DEFAULT=-1, PETSCDTNODES_GAUSSJACOBI, PETSCDTNODES_EQUISPACED, PETSCDTNODES_TANHSINH} PetscDTNodeType;

 49: PETSC_EXTERN const char *const PetscDTNodeTypes[];

 51: PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
 52: PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
 53: PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
 54: PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
 55: PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
 56: PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
 57: PETSC_EXTERN PetscErrorCode PetscQuadratureEqual(PetscQuadrature, PetscQuadrature, PetscBool*);
 58: PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
 59: PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
 60: PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
 61: PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);

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

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

 68: PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
 69: PETSC_EXTERN PetscErrorCode PetscDTJacobiNorm(PetscReal,PetscReal,PetscInt,PetscReal *);
 70: PETSC_EXTERN PetscErrorCode PetscDTJacobiEval(PetscInt,PetscReal,PetscReal,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
 71: PETSC_EXTERN PetscErrorCode PetscDTJacobiEvalJet(PetscReal,PetscReal,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]);
 72: PETSC_EXTERN PetscErrorCode PetscDTPKDEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscReal[]);
 73: PETSC_EXTERN PetscErrorCode PetscDTPTrimmedSize(PetscInt,PetscInt,PetscInt,PetscInt*);
 74: PETSC_EXTERN PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt,PetscInt,const PetscReal[],PetscInt,PetscInt,PetscInt,PetscReal[]);
 75: PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
 76: PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
 77: PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*);
 78: PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
 79: PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
 80: PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
 81: PETSC_EXTERN PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);

 83: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
 84: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);
 85: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(const PetscReal[], void *, PetscReal *), PetscReal, PetscReal, PetscInt, void *, PetscReal *);

 87: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
 88: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 89: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 90: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
 91: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
 92: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 93: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 94: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 95: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);

 97: PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
 98: PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
 99: PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
100: PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
101: PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
102: PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
103: PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
104: PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
105: PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);

107: PETSC_EXTERN PetscErrorCode PetscDTBaryToIndex(PetscInt,PetscInt,const PetscInt[],PetscInt*);
108: PETSC_EXTERN PetscErrorCode PetscDTIndexToBary(PetscInt,PetscInt,PetscInt,PetscInt[]);
109: PETSC_EXTERN PetscErrorCode PetscDTGradedOrderToIndex(PetscInt,const PetscInt[],PetscInt*);
110: PETSC_EXTERN PetscErrorCode PetscDTIndexToGradedOrder(PetscInt,PetscInt,PetscInt[]);

112: #if defined(PETSC_USE_64BIT_INDICES)
113: #define PETSC_FACTORIAL_MAX 20
114: #define PETSC_BINOMIAL_MAX  61
115: #else
116: #define PETSC_FACTORIAL_MAX 12
117: #define PETSC_BINOMIAL_MAX  29
118: #endif

120: /*MC
121:    PetscDTFactorial - Approximate n! as a real number

123:    Input Parameter:
124: .  n - a non-negative integer

126:    Output Parameter:
127: .  factorial - n!

129:    Level: beginner
130: M*/
131: static inline PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
132: {
133:   PetscReal f = 1.0;

135:   *factorial = -1.0;
137:   for (PetscInt i = 1; i < n+1; ++i) f *= (PetscReal)i;
138:   *factorial = f;
139:   return 0;
140: }

142: /*MC
143:    PetscDTFactorialInt - Compute n! as an integer

145:    Input Parameter:
146: .  n - a non-negative integer

148:    Output Parameter:
149: .  factorial - n!

151:    Level: beginner

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

159:   *factorial = -1;
161:   if (n <= 12) {
162:     *factorial = facLookup[n];
163:   } else {
164:     PetscInt f = facLookup[12];
165:     PetscInt i;

167:     for (i = 13; i < n+1; ++i) f *= i;
168:     *factorial = f;
169:   }
170:   return 0;
171: }

173: /*MC
174:    PetscDTBinomial - Approximate the binomial coefficient "n choose k"

176:    Input Parameters:
177: +  n - a non-negative integer
178: -  k - an integer between 0 and n, inclusive

180:    Output Parameter:
181: .  binomial - approximation of the binomial coefficient n choose k

183:    Level: beginner
184: M*/
185: static inline PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
186: {
188:   *binomial = -1.0;
190:   if (n <= 3) {
191:     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};

193:     *binomial = (PetscReal)binomLookup[n][k];
194:   } else {
195:     PetscReal binom = 1.0;

197:     k = PetscMin(k, n - k);
198:     for (PetscInt i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
199:     *binomial = binom;
200:   }
201:   return 0;
202: }

204: /*MC
205:    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"

207:    Input Parameters:
208: +  n - a non-negative integer
209: -  k - an integer between 0 and n, inclusive

211:    Output Parameter:
212: .  binomial - the binomial coefficient n choose k

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

216:    Level: beginner
217: M*/
218: static inline PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
219: {
220:   PetscInt bin;

222:   *binomial = -1;
225:   if (n <= 3) {
226:     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};

228:     bin = binomLookup[n][k];
229:   } else {
230:     PetscInt binom = 1;

232:     k = PetscMin(k, n - k);
233:     for (PetscInt i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
234:     bin = binom;
235:   }
236:   *binomial = bin;
237:   return 0;
238: }

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

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

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

253:    Output Parameters:
254: +  perm - the permuted list of the integers [0, ..., n-1]
255: -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.

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

259:    Level: beginner
260: M*/
261: static inline PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
262: {
263:   PetscInt  odd = 0;
264:   PetscInt  i;
265:   PetscInt  work[PETSC_FACTORIAL_MAX];
266:   PetscInt *w;

268:   if (isOdd) *isOdd = PETSC_FALSE;
270:   w = &work[n - 2];
271:   for (i = 2; i <= n; i++) {
272:     *(w--) = k % i;
273:     k /= i;
274:   }
275:   for (i = 0; i < n; i++) perm[i] = i;
276:   for (i = 0; i < n - 1; i++) {
277:     PetscInt s = work[i];
278:     PetscInt swap = perm[i];

280:     perm[i] = perm[i + s];
281:     perm[i + s] = swap;
282:     odd ^= (!!s);
283:   }
284:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
285:   return 0;
286: }

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

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

295:    Output Parameters:
296: +  k - an integer in [0, n!)
297: -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.

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

301:    Level: beginner
302: M*/
303: static inline PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
304: {
305:   PetscInt  odd = 0;
306:   PetscInt  i, idx;
307:   PetscInt  work[PETSC_FACTORIAL_MAX];
308:   PetscInt  iwork[PETSC_FACTORIAL_MAX];

311:   *k = -1;
312:   if (isOdd) *isOdd = PETSC_FALSE;
314:   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
315:   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
316:   for (idx = 0, i = 0; i < n - 1; i++) {
317:     PetscInt j = perm[i];
318:     PetscInt icur = work[i];
319:     PetscInt jloc = iwork[j];
320:     PetscInt diff = jloc - i;

322:     idx = idx * (n - i) + diff;
323:     /* swap (i, jloc) */
324:     work[i] = j;
325:     work[jloc] = icur;
326:     iwork[j] = i;
327:     iwork[icur] = jloc;
328:     odd ^= (!!diff);
329:   }
330:   *k = idx;
331:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
332:   return 0;
333: }

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

339:    Input Parameters:
340: +  n - a non-negative integer (see note about limits below)
341: .  k - an integer in [0, n]
342: -  j - an index in [0, n choose k)

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

347:    Note: this is limited by arguments such that n choose k can be represented by PetscInt

349:    Level: beginner

351: .seealso: PetscDTSubsetIndex()
352: M*/
353: static inline PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
354: {
355:   PetscInt Nk;

358:   PetscDTBinomialInt(n, k, &Nk);
359:   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
360:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
361:     PetscInt Nminusk = Nk - Nminuskminus;

363:     if (j < Nminuskminus) {
364:       subset[l++] = i;
365:       Nk = Nminuskminus;
366:     } else {
367:       j -= Nminuskminus;
368:       Nk = Nminusk;
369:     }
370:   }
371:   return 0;
372: }

374: /*MC
375:    PetscDTSubsetIndex - Convert an ordered subset of k integers from the set [0, ..., n - 1] to its encoding as an integers in [0, n choose k) in lexicographic order.  This is the inverse of PetscDTEnumSubset.

377:    Input Parameters:
378: +  n - a non-negative integer (see note about limits below)
379: .  k - an integer in [0, n]
380: -  subset - an ordered subset of the integers [0, ..., n - 1]

382:    Output Parameter:
383: .  index - the rank of the subset in lexicographic order

385:    Note: this is limited by arguments such that n choose k can be represented by PetscInt

387:    Level: beginner

389: .seealso: PetscDTEnumSubset()
390: M*/
391: static inline PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
392: {
393:   PetscInt j = 0, Nk;

395:   *index = -1;
396:   PetscDTBinomialInt(n, k, &Nk);
397:   for (PetscInt i = 0, l = 0; i < n && l < k; i++) {
398:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
399:     PetscInt Nminusk = Nk - Nminuskminus;

401:     if (subset[l] == i) {
402:       l++;
403:       Nk = Nminuskminus;
404:     } else {
405:       j += Nminuskminus;
406:       Nk = Nminusk;
407:     }
408:   }
409:   *index = j;
410:   return 0;
411: }

413: /*MC
414:    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.

416:    Input Parameters:
417: +  n - a non-negative integer (see note about limits below)
418: .  k - an integer in [0, n]
419: -  j - an index in [0, n choose k)

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

425:    Note: this is limited by arguments such that n choose k can be represented by PetscInt

427:    Level: beginner

429: .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
430: M*/
431: static inline PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
432: {
433:   PetscInt i, l, m, Nk, odd = 0;
434:   PetscInt *subcomp = perm+k;

436:   if (isOdd) *isOdd = PETSC_FALSE;
437:   PetscDTBinomialInt(n, k, &Nk);
438:   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
439:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
440:     PetscInt Nminusk = Nk - Nminuskminus;

442:     if (j < Nminuskminus) {
443:       perm[l++] = i;
444:       Nk = Nminuskminus;
445:     } else {
446:       subcomp[m++] = i;
447:       j -= Nminuskminus;
448:       odd ^= ((k - l) & 1);
449:       Nk = Nminusk;
450:     }
451:   }
452:   for (; i < n; i++) subcomp[m++] = i;
453:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
454:   return 0;
455: }

457: struct _p_PetscTabulation {
458:   PetscInt    K;    /* Indicates a k-jet, namely tabulated derivatives up to order k */
459:   PetscInt    Nr;   /* The number of tabulation replicas (often 1) */
460:   PetscInt    Np;   /* The number of tabulation points in a replica */
461:   PetscInt    Nb;   /* The number of functions tabulated */
462:   PetscInt    Nc;   /* The number of function components */
463:   PetscInt    cdim; /* The coordinate dimension */
464:   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
465:                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
466:                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
467:                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
468: };
469: typedef struct _p_PetscTabulation *PetscTabulation;

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

473: typedef enum {DTPROB_DENSITY_CONSTANT, DTPROB_DENSITY_GAUSSIAN, DTPROB_DENSITY_MAXWELL_BOLTZMANN, DTPROB_NUM_DENSITY} DTProbDensityType;
474: PETSC_EXTERN const char * const DTProbDensityTypes[];

476: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
477: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal[], const PetscReal[], PetscReal[]);
478: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
479: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal[], const PetscReal[], PetscReal[]);
480: PETSC_EXTERN PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
481: PETSC_EXTERN PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal[], const PetscReal[], PetscReal[]);
482: PETSC_EXTERN PetscErrorCode PetscPDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
483: PETSC_EXTERN PetscErrorCode PetscCDFGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
484: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal[], const PetscReal[], PetscReal[]);
485: PETSC_EXTERN PetscErrorCode PetscPDFGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
486: PETSC_EXTERN PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal[], const PetscReal[], PetscReal[]);
487: PETSC_EXTERN PetscErrorCode PetscPDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
488: PETSC_EXTERN PetscErrorCode PetscCDFConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
489: PETSC_EXTERN PetscErrorCode PetscPDFSampleConstant1D(const PetscReal[], const PetscReal[], PetscReal[]);
490: PETSC_EXTERN PetscErrorCode PetscProbCreateFromOptions(PetscInt, const char[], const char[], PetscProbFunc *, PetscProbFunc *, PetscProbFunc *);

492: #include <petscvec.h>

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

496: #endif