Actual source code: petscdt.h

petsc-master 2020-01-17
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.;
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 *= 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:   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);
157:   if (n <= 3) {
158:     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};

160:     *binomial = binomLookup[n][k];
161:   } else {
162:     PetscReal binom = 1.;
163:     PetscInt  i;

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

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

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

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

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

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

191:   *binomial = -1;
192:   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);
193:   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);
194:   if (n <= 3) {
195:     PetscInt binomLookup[4][4] = {{1, 0, 0, 0}, {1, 1, 0, 0}, {1, 2, 1, 0}, {1, 3, 3, 1}};

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

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

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

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

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

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

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

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

239:   if (isOdd) *isOdd = PETSC_FALSE;
240:   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);
241:   w = &work[n - 2];
242:   for (i = 2; i <= n; i++) {
243:     *(w--) = k % i;
244:     k /= i;
245:   }
246:   for (i = 0; i < n; i++) perm[i] = i;
247:   for (i = 0; i < n - 1; i++) {
248:     PetscInt s = work[i];
249:     PetscInt swap = perm[i];

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

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

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

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

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

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

282:   *k = -1;
283:   if (isOdd) *isOdd = PETSC_FALSE;
284:   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);
285:   for (i = 0; i < n; i++) work[i] = i;  /* partial permutation */
286:   for (i = 0; i < n; i++) iwork[i] = i; /* partial permutation inverse */
287:   for (idx = 0, i = 0; i < n - 1; i++) {
288:     PetscInt j = perm[i];
289:     PetscInt icur = work[i];
290:     PetscInt jloc = iwork[j];
291:     PetscInt diff = jloc - i;

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

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

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

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

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

320:    Level: beginner

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

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

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

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

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

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

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

359:    Level: beginner

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

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

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

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

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

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

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

401:    Level: beginner

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

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

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

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

451: #endif