Actual source code: petscdt.h

petsc-master 2020-02-18
Report Typos and Errors
  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: PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
 32: PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
 33: PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
 34: PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
 35: PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
 36: PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
 37: PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
 38: PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
 39: PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
 40: PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);

 42: PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);

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

 46: PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
 47: PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
 48: PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGaussLobattoLegendreCreateType,PetscReal*,PetscReal*);
 49: PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
 50: PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
 51: PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);

 53: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
 54: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
 55: PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);

 57: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt, PetscReal *, PetscReal *, const PetscReal *, PetscReal *);
 58: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 59: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 60: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
 61: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***, PetscReal ***);
 62: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 63: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 64: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt, PetscReal *, PetscReal *, PetscReal ***);
 65: PETSC_EXTERN PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt, PetscReal *, PetscReal *, PetscReal ***);

 67: PETSC_EXTERN PetscErrorCode PetscDTAltVApply(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
 68: PETSC_EXTERN PetscErrorCode PetscDTAltVWedge(PetscInt, PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
 69: PETSC_EXTERN PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);
 70: PETSC_EXTERN PetscErrorCode PetscDTAltVPullback(PetscInt, PetscInt, const PetscReal *, PetscInt, const PetscReal *, PetscReal *);
 71: PETSC_EXTERN PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt, PetscInt, const PetscReal *, PetscInt, PetscReal *);
 72: PETSC_EXTERN PetscErrorCode PetscDTAltVInterior(PetscInt, PetscInt, const PetscReal *, const PetscReal *, PetscReal *);
 73: PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt, PetscInt, const PetscReal *, PetscReal *);
 74: PETSC_EXTERN PetscErrorCode PetscDTAltVInteriorPattern(PetscInt, PetscInt, PetscInt (*)[3]);
 75: PETSC_EXTERN PetscErrorCode PetscDTAltVStar(PetscInt, PetscInt, PetscInt, const PetscReal *, PetscReal *);

 77: #if defined(PETSC_USE_64BIT_INDICES)
 78: #define PETSC_FACTORIAL_MAX 20
 79: #define PETSC_BINOMIAL_MAX  61
 80: #else
 81: #define PETSC_FACTORIAL_MAX 12
 82: #define PETSC_BINOMIAL_MAX  29
 83: #endif

 85: /*MC
 86:    PetscDTFactorial - Approximate n! as a real number

 88:    Input Arguments:
 89: .  n - a non-negative integer

 91:    Output Arguments:
 92: .  factorial - n!

 94:    Level: beginner
 95: M*/
 96: PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial(PetscInt n, PetscReal *factorial)
 97: {
 98:   PetscReal f = 1.0;
 99:   PetscInt  i;

102:   *factorial = -1.0;
103:   if (n < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Factorial called with negative number %D\n", n);
104:   for (i = 1; i < n+1; ++i) f *= (PetscReal)i;
105:   *factorial = f;
106:   return(0);
107: }

109: /*MC
110:    PetscDTFactorialInt - Compute n! as an integer

112:    Input Arguments:
113: .  n - a non-negative integer

115:    Output Arguments:
116: .  factorial - n!

118:    Level: beginner

120:    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.
121: M*/
122: PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorialInt(PetscInt n, PetscInt *factorial)
123: {
124:   PetscInt facLookup[13] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600};

127:   *factorial = -1;
128:   if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX);
129:   if (n <= 12) {
130:     *factorial = facLookup[n];
131:   } else {
132:     PetscInt f = facLookup[12];
133:     PetscInt i;

135:     for (i = 13; i < n+1; ++i) f *= i;
136:     *factorial = f;
137:   }
138:   return(0);
139: }

141: /*MC
142:    PetscDTBinomial - Approximate the binomial coefficient "n choose k"

144:    Input Arguments:
145: +  n - a non-negative integer
146: -  k - an integer between 0 and n, inclusive

148:    Output Arguments:
149: .  binomial - approximation of the binomial coefficient n choose k

151:    Level: beginner
152: M*/
153: PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomial(PetscInt n, PetscInt k, PetscReal *binomial)
154: {
156:   *binomial = -1.0;
157:   if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n\n", n, k);
158:   if (n <= 3) {
159:     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};

161:     *binomial = (PetscReal)binomLookup[n][k];
162:   } else {
163:     PetscReal binom = 1.0;
164:     PetscInt  i;

166:     k = PetscMin(k, n - k);
167:     for (i = 0; i < k; i++) binom = (binom * (PetscReal)(n - i)) / (PetscReal)(i + 1);
168:     *binomial = binom;
169:   }
170:   return(0);
171: }

173: /*MC
174:    PetscDTBinomialInt - Compute the binomial coefficient "n choose k"

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

180:    Output Arguments:
181: .  binomial - the binomial coefficient n choose k

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

185:    Level: beginner
186: M*/
187: PETSC_STATIC_INLINE PetscErrorCode PetscDTBinomialInt(PetscInt n, PetscInt k, PetscInt *binomial)
188: {
189:   PetscInt bin;

192:   *binomial = -1;
193:   if (n < 0 || k < 0 || k > n) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial arguments (%D %D) must be non-negative, k <= n\n", n, k);
194:   if (n > PETSC_BINOMIAL_MAX) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Binomial elements %D is larger than max for PetscInt, %D\n", n, PETSC_BINOMIAL_MAX);
195:   if (n <= 3) {
196:     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};

198:     bin = binomLookup[n][k];
199:   } else {
200:     PetscInt  binom = 1;
201:     PetscInt  i;

203:     k = PetscMin(k, n - k);
204:     for (i = 0; i < k; i++) binom = (binom * (n - i)) / (i + 1);
205:     bin = binom;
206:   }
207:   *binomial = bin;
208:   return(0);
209: }

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

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

220:    Input Arguments:
221: +  n - a non-negative integer (see note about limits below)
222: -  k - an integer in [0, n!)

224:    Output Arguments:
225: +  perm - the permuted list of the integers [0, ..., n-1]
226: -  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.

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

230:    Level: beginner
231: M*/
232: PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumPerm(PetscInt n, PetscInt k, PetscInt *perm, PetscBool *isOdd)
233: {
234:   PetscInt  odd = 0;
235:   PetscInt  i;
236:   PetscInt  work[PETSC_FACTORIAL_MAX];
237:   PetscInt *w;

240:   if (isOdd) *isOdd = PETSC_FALSE;
241:   if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX);
242:   w = &work[n - 2];
243:   for (i = 2; i <= n; i++) {
244:     *(w--) = k % i;
245:     k /= i;
246:   }
247:   for (i = 0; i < n; i++) perm[i] = i;
248:   for (i = 0; i < n - 1; i++) {
249:     PetscInt s = work[i];
250:     PetscInt swap = perm[i];

252:     perm[i] = perm[i + s];
253:     perm[i + s] = swap;
254:     odd ^= (!!s);
255:   }
256:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
257:   return(0);
258: }

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

263:    Input Arguments:
264: +  n - a non-negative integer (see note about limits below)
265: -  perm - the permuted list of the integers [0, ..., n-1]

267:    Output Arguments:
268: +  k - an integer in [0, n!)
269: .  isOdd - if not NULL, returns wether the permutation used an even or odd number of swaps.

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

273:    Level: beginner
274: M*/
275: PETSC_STATIC_INLINE PetscErrorCode PetscDTPermIndex(PetscInt n, const PetscInt *perm, PetscInt *k, PetscBool *isOdd)
276: {
277:   PetscInt  odd = 0;
278:   PetscInt  i, idx;
279:   PetscInt  work[PETSC_FACTORIAL_MAX];
280:   PetscInt  iwork[PETSC_FACTORIAL_MAX];

283:   *k = -1;
284:   if (isOdd) *isOdd = PETSC_FALSE;
285:   if (n < 0 || n > PETSC_FACTORIAL_MAX) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of elements %D is not in supported range [0,%D]\n",n,PETSC_FACTORIAL_MAX);
286:   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
287:   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
288:   for (idx = 0, i = 0; i < n - 1; i++) {
289:     PetscInt j = perm[i];
290:     PetscInt icur = work[i];
291:     PetscInt jloc = iwork[j];
292:     PetscInt diff = jloc - i;

294:     idx = idx * (n - i) + diff;
295:     /* swap (i, jloc) */
296:     work[i] = j;
297:     work[jloc] = icur;
298:     iwork[j] = i;
299:     iwork[icur] = jloc;
300:     odd ^= (!!diff);
301:   }
302:   *k = idx;
303:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
304:   return(0);
305: }

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

311:    Input Arguments:
312: +  n - a non-negative integer (see note about limits below)
313: .  k - an integer in [0, n]
314: -  j - an index in [0, n choose k)

316:    Output Arguments:
317: .  subset - the jth subset of size k of the integers [0, ..., n - 1]

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

321:    Level: beginner

323: .seealso: PetscDTSubsetIndex()
324: M*/
325: PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSubset(PetscInt n, PetscInt k, PetscInt j, PetscInt *subset)
326: {
327:   PetscInt       Nk, i, l;

331:   PetscDTBinomialInt(n, k, &Nk);
332:   for (i = 0, l = 0; i < n && l < k; i++) {
333:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
334:     PetscInt Nminusk = Nk - Nminuskminus;

336:     if (j < Nminuskminus) {
337:       subset[l++] = i;
338:       Nk = Nminuskminus;
339:     } else {
340:       j -= Nminuskminus;
341:       Nk = Nminusk;
342:     }
343:   }
344:   return(0);
345: }

347: /*MC
348:    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.

350:    Input Arguments:
351: +  n - a non-negative integer (see note about limits below)
352: .  k - an integer in [0, n]
353: -  subset - an ordered subset of the integers [0, ..., n - 1]

355:    Output Arguments:
356: .  index - the rank of the subset in lexicographic order

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

360:    Level: beginner

362: .seealso: PetscDTEnumSubset()
363: M*/
364: PETSC_STATIC_INLINE PetscErrorCode PetscDTSubsetIndex(PetscInt n, PetscInt k, const PetscInt *subset, PetscInt *index)
365: {
366:   PetscInt       i, j = 0, l, Nk;

370:   *index = -1;
371:   PetscDTBinomialInt(n, k, &Nk);
372:   for (i = 0, l = 0; i < n && l < k; i++) {
373:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
374:     PetscInt Nminusk = Nk - Nminuskminus;

376:     if (subset[l] == i) {
377:       l++;
378:       Nk = Nminuskminus;
379:     } else {
380:       j += Nminuskminus;
381:       Nk = Nminusk;
382:     }
383:   }
384:   *index = j;
385:   return(0);
386: }

388: /*MC
389:    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.

391:    Input Arguments:
392: +  n - a non-negative integer (see note about limits below)
393: .  k - an integer in [0, n]
394: -  j - an index in [0, n choose k)

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

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

402:    Level: beginner

404: .seealso: PetscDTEnumSubset(), PetscDTSubsetIndex()
405: M*/
406: PETSC_STATIC_INLINE PetscErrorCode PetscDTEnumSplit(PetscInt n, PetscInt k, PetscInt j, PetscInt *perm, PetscBool *isOdd)
407: {
408:   PetscInt       i, l, m, *subcomp, Nk;
409:   PetscInt       odd;

413:   if (isOdd) *isOdd = PETSC_FALSE;
414:   PetscDTBinomialInt(n, k, &Nk);
415:   odd = 0;
416:   subcomp = &perm[k];
417:   for (i = 0, l = 0, m = 0; i < n && l < k; i++) {
418:     PetscInt Nminuskminus = (Nk * (k - l)) / (n - i);
419:     PetscInt Nminusk = Nk - Nminuskminus;

421:     if (j < Nminuskminus) {
422:       perm[l++] = i;
423:       Nk = Nminuskminus;
424:     } else {
425:       subcomp[m++] = i;
426:       j -= Nminuskminus;
427:       odd ^= ((k - l) & 1);
428:       Nk = Nminusk;
429:     }
430:   }
431:   for (; i < n; i++) {
432:     subcomp[m++] = i;
433:   }
434:   if (isOdd) *isOdd = odd ? PETSC_TRUE : PETSC_FALSE;
435:   return(0);
436: }

438: struct _p_PetscTabulation {
439:   PetscInt    K;    /* Indicates a k-jet, namely tabulated derviatives up to order k */
440:   PetscInt    Nr;   /* THe number of tabulation replicas (often 1) */
441:   PetscInt    Np;   /* The number of tabulation points in a replica */
442:   PetscInt    Nb;   /* The number of functions tabulated */
443:   PetscInt    Nc;   /* The number of function components */
444:   PetscInt    cdim; /* The coordinate dimension */
445:   PetscReal **T;    /* The tabulation T[K] of functions and their derivatives
446:                        T[0] = B[Nr*Np][Nb][Nc]:             The basis function values at quadrature points
447:                        T[1] = D[Nr*Np][Nb][Nc][cdim]:       The basis function derivatives at quadrature points
448:                        T[2] = H[Nr*Np][Nb][Nc][cdim][cdim]: The basis function second derivatives at quadrature points */
449: };
450: typedef struct _p_PetscTabulation *PetscTabulation;

452: #endif