Actual source code: dt.c

  1: /* Discretization tools */

  3: #include <petscdt.h>
  4: #include <petscblaslapack.h>
  5: #include <petsc/private/petscimpl.h>
  6: #include <petsc/private/dtimpl.h>
  7: #include <petsc/private/petscfeimpl.h>
  8: #include <petscviewer.h>
  9: #include <petscdmplex.h>
 10: #include <petscdmshell.h>

 12: #if defined(PETSC_HAVE_MPFR)
 13:   #include <mpfr.h>
 14: #endif

 16: const char *const        PetscDTNodeTypes_shifted[] = {"default", "gaussjacobi", "equispaced", "tanhsinh", "PetscDTNodeType", "PETSCDTNODES_", NULL};
 17: const char *const *const PetscDTNodeTypes           = PetscDTNodeTypes_shifted + 1;

 19: const char *const        PetscDTSimplexQuadratureTypes_shifted[] = {"default", "conic", "minsym", "diagsym", "PetscDTSimplexQuadratureType", "PETSCDTSIMPLEXQUAD_", NULL};
 20: const char *const *const PetscDTSimplexQuadratureTypes           = PetscDTSimplexQuadratureTypes_shifted + 1;

 22: static PetscBool GolubWelschCite       = PETSC_FALSE;
 23: const char       GolubWelschCitation[] = "@article{GolubWelsch1969,\n"
 24:                                          "  author  = {Golub and Welsch},\n"
 25:                                          "  title   = {Calculation of Quadrature Rules},\n"
 26:                                          "  journal = {Math. Comp.},\n"
 27:                                          "  volume  = {23},\n"
 28:                                          "  number  = {106},\n"
 29:                                          "  pages   = {221--230},\n"
 30:                                          "  year    = {1969}\n}\n";

 32: /* Numerical tests in src/dm/dt/tests/ex1.c show that when computing the nodes and weights of Gauss-Jacobi
 33:    quadrature rules:

 35:    - in double precision, Newton's method and Golub & Welsch both work for moderate degrees (< 100),
 36:    - in single precision, Newton's method starts producing incorrect roots around n = 15, but
 37:      the weights from Golub & Welsch become a problem before then: they produces errors
 38:      in computing the Jacobi-polynomial Gram matrix around n = 6.

 40:    So we default to Newton's method (required fewer dependencies) */
 41: PetscBool PetscDTGaussQuadratureNewton_Internal = PETSC_TRUE;

 43: PetscClassId PETSCQUADRATURE_CLASSID = 0;

 45: /*@
 46:   PetscQuadratureCreate - Create a `PetscQuadrature` object

 48:   Collective

 50:   Input Parameter:
 51: . comm - The communicator for the `PetscQuadrature` object

 53:   Output Parameter:
 54: . q - The `PetscQuadrature` object

 56:   Level: beginner

 58: .seealso: `PetscQuadrature`, `Petscquadraturedestroy()`, `PetscQuadratureGetData()`
 59: @*/
 60: PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
 61: {
 62:   PetscFunctionBegin;
 63:   PetscAssertPointer(q, 2);
 64:   PetscCall(DMInitializePackage());
 65:   PetscCall(PetscHeaderCreate(*q, PETSCQUADRATURE_CLASSID, "PetscQuadrature", "Quadrature", "DT", comm, PetscQuadratureDestroy, PetscQuadratureView));
 66:   (*q)->ct        = DM_POLYTOPE_UNKNOWN;
 67:   (*q)->dim       = -1;
 68:   (*q)->Nc        = 1;
 69:   (*q)->order     = -1;
 70:   (*q)->numPoints = 0;
 71:   (*q)->points    = NULL;
 72:   (*q)->weights   = NULL;
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }

 76: /*@
 77:   PetscQuadratureDuplicate - Create a deep copy of the `PetscQuadrature` object

 79:   Collective

 81:   Input Parameter:
 82: . q - The `PetscQuadrature` object

 84:   Output Parameter:
 85: . r - The new `PetscQuadrature` object

 87:   Level: beginner

 89: .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`, `PetscQuadratureDestroy()`, `PetscQuadratureGetData()`
 90: @*/
 91: PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
 92: {
 93:   DMPolytopeType   ct;
 94:   PetscInt         order, dim, Nc, Nq;
 95:   const PetscReal *points, *weights;
 96:   PetscReal       *p, *w;

 98:   PetscFunctionBegin;
 99:   PetscAssertPointer(q, 1);
100:   PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)q), r));
101:   PetscCall(PetscQuadratureGetCellType(q, &ct));
102:   PetscCall(PetscQuadratureSetCellType(*r, ct));
103:   PetscCall(PetscQuadratureGetOrder(q, &order));
104:   PetscCall(PetscQuadratureSetOrder(*r, order));
105:   PetscCall(PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights));
106:   PetscCall(PetscMalloc1(Nq * dim, &p));
107:   PetscCall(PetscMalloc1(Nq * Nc, &w));
108:   PetscCall(PetscArraycpy(p, points, Nq * dim));
109:   PetscCall(PetscArraycpy(w, weights, Nc * Nq));
110:   PetscCall(PetscQuadratureSetData(*r, dim, Nc, Nq, p, w));
111:   PetscFunctionReturn(PETSC_SUCCESS);
112: }

114: /*@
115:   PetscQuadratureDestroy - Destroys a `PetscQuadrature` object

117:   Collective

119:   Input Parameter:
120: . q - The `PetscQuadrature` object

122:   Level: beginner

124: .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`, `PetscQuadratureGetData()`
125: @*/
126: PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
127: {
128:   PetscFunctionBegin;
129:   if (!*q) PetscFunctionReturn(PETSC_SUCCESS);
131:   if (--((PetscObject)*q)->refct > 0) {
132:     *q = NULL;
133:     PetscFunctionReturn(PETSC_SUCCESS);
134:   }
135:   PetscCall(PetscFree((*q)->points));
136:   PetscCall(PetscFree((*q)->weights));
137:   PetscCall(PetscHeaderDestroy(q));
138:   PetscFunctionReturn(PETSC_SUCCESS);
139: }

141: /*@
142:   PetscQuadratureGetCellType - Return the cell type of the integration domain

144:   Not Collective

146:   Input Parameter:
147: . q - The `PetscQuadrature` object

149:   Output Parameter:
150: . ct - The cell type of the integration domain

152:   Level: intermediate

154: .seealso: `PetscQuadrature`, `PetscQuadratureSetCellType()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
155: @*/
156: PetscErrorCode PetscQuadratureGetCellType(PetscQuadrature q, DMPolytopeType *ct)
157: {
158:   PetscFunctionBegin;
160:   PetscAssertPointer(ct, 2);
161:   *ct = q->ct;
162:   PetscFunctionReturn(PETSC_SUCCESS);
163: }

165: /*@
166:   PetscQuadratureSetCellType - Set the cell type of the integration domain

168:   Not Collective

170:   Input Parameters:
171: + q  - The `PetscQuadrature` object
172: - ct - The cell type of the integration domain

174:   Level: intermediate

176: .seealso: `PetscQuadrature`, `PetscQuadratureGetCellType()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
177: @*/
178: PetscErrorCode PetscQuadratureSetCellType(PetscQuadrature q, DMPolytopeType ct)
179: {
180:   PetscFunctionBegin;
182:   q->ct = ct;
183:   PetscFunctionReturn(PETSC_SUCCESS);
184: }

186: /*@
187:   PetscQuadratureGetOrder - Return the order of the method in the `PetscQuadrature`

189:   Not Collective

191:   Input Parameter:
192: . q - The `PetscQuadrature` object

194:   Output Parameter:
195: . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated

197:   Level: intermediate

199: .seealso: `PetscQuadrature`, `PetscQuadratureSetOrder()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
200: @*/
201: PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
202: {
203:   PetscFunctionBegin;
205:   PetscAssertPointer(order, 2);
206:   *order = q->order;
207:   PetscFunctionReturn(PETSC_SUCCESS);
208: }

210: /*@
211:   PetscQuadratureSetOrder - Set the order of the method in the `PetscQuadrature`

213:   Not Collective

215:   Input Parameters:
216: + q     - The `PetscQuadrature` object
217: - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated

219:   Level: intermediate

221: .seealso: `PetscQuadrature`, `PetscQuadratureGetOrder()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
222: @*/
223: PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
224: {
225:   PetscFunctionBegin;
227:   q->order = order;
228:   PetscFunctionReturn(PETSC_SUCCESS);
229: }

231: /*@
232:   PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated

234:   Not Collective

236:   Input Parameter:
237: . q - The `PetscQuadrature` object

239:   Output Parameter:
240: . Nc - The number of components

242:   Level: intermediate

244:   Note:
245:   We are performing an integral $\int f(x) w(x) dx$, where both $f$ and $w$ (the weight) have `Nc` components.

247: .seealso: `PetscQuadrature`, `PetscQuadratureSetNumComponents()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
248: @*/
249: PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc)
250: {
251:   PetscFunctionBegin;
253:   PetscAssertPointer(Nc, 2);
254:   *Nc = q->Nc;
255:   PetscFunctionReturn(PETSC_SUCCESS);
256: }

258: /*@
259:   PetscQuadratureSetNumComponents - Sets the number of components for functions to be integrated

261:   Not Collective

263:   Input Parameters:
264: + q  - The `PetscQuadrature` object
265: - Nc - The number of components

267:   Level: intermediate

269:   Note:
270:   We are performing an integral $\int f(x) w(x) dx$, where both $f$ and $w$ (the weight) have `Nc` components.

272: .seealso: `PetscQuadrature`, `PetscQuadratureGetNumComponents()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
273: @*/
274: PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc)
275: {
276:   PetscFunctionBegin;
278:   q->Nc = Nc;
279:   PetscFunctionReturn(PETSC_SUCCESS);
280: }

282: /*@C
283:   PetscQuadratureGetData - Returns the data defining the `PetscQuadrature`

285:   Not Collective

287:   Input Parameter:
288: . q - The `PetscQuadrature` object

290:   Output Parameters:
291: + dim     - The spatial dimension
292: . Nc      - The number of components
293: . npoints - The number of quadrature points
294: . points  - The coordinates of each quadrature point
295: - weights - The weight of each quadrature point

297:   Level: intermediate

299:   Note:
300:   All output arguments are optional, pass `NULL` for any argument not required

302:   Fortran Note:
303:   Call `PetscQuadratureRestoreData()` when you are done with the data

305: .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`, `PetscQuadratureSetData()`
306: @*/
307: PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PeOp PetscInt *dim, PeOp PetscInt *Nc, PeOp PetscInt *npoints, PeOp const PetscReal *points[], PeOp const PetscReal *weights[])
308: {
309:   PetscFunctionBegin;
311:   if (dim) {
312:     PetscAssertPointer(dim, 2);
313:     *dim = q->dim;
314:   }
315:   if (Nc) {
316:     PetscAssertPointer(Nc, 3);
317:     *Nc = q->Nc;
318:   }
319:   if (npoints) {
320:     PetscAssertPointer(npoints, 4);
321:     *npoints = q->numPoints;
322:   }
323:   if (points) {
324:     PetscAssertPointer(points, 5);
325:     *points = q->points;
326:   }
327:   if (weights) {
328:     PetscAssertPointer(weights, 6);
329:     *weights = q->weights;
330:   }
331:   PetscFunctionReturn(PETSC_SUCCESS);
332: }

334: /*@
335:   PetscQuadratureEqual - determine whether two quadratures are equivalent

337:   Input Parameters:
338: + A - A `PetscQuadrature` object
339: - B - Another `PetscQuadrature` object

341:   Output Parameter:
342: . equal - `PETSC_TRUE` if the quadratures are the same

344:   Level: intermediate

346: .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`
347: @*/
348: PetscErrorCode PetscQuadratureEqual(PetscQuadrature A, PetscQuadrature B, PetscBool *equal)
349: {
350:   PetscFunctionBegin;
353:   PetscAssertPointer(equal, 3);
354:   *equal = PETSC_FALSE;
355:   if (A->ct != B->ct || A->dim != B->dim || A->Nc != B->Nc || A->order != B->order || A->numPoints != B->numPoints) PetscFunctionReturn(PETSC_SUCCESS);
356:   for (PetscInt i = 0; i < A->numPoints * A->dim; i++) {
357:     if (PetscAbsReal(A->points[i] - B->points[i]) > PETSC_SMALL) PetscFunctionReturn(PETSC_SUCCESS);
358:   }
359:   if (!A->weights && !B->weights) {
360:     *equal = PETSC_TRUE;
361:     PetscFunctionReturn(PETSC_SUCCESS);
362:   }
363:   if (A->weights && B->weights) {
364:     for (PetscInt i = 0; i < A->numPoints; i++) {
365:       if (PetscAbsReal(A->weights[i] - B->weights[i]) > PETSC_SMALL) PetscFunctionReturn(PETSC_SUCCESS);
366:     }
367:     *equal = PETSC_TRUE;
368:   }
369:   PetscFunctionReturn(PETSC_SUCCESS);
370: }

372: static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[])
373: {
374:   PetscScalar *Js, *Jinvs;
375:   PetscInt     i, j, k;
376:   PetscBLASInt bm, bn, info;

378:   PetscFunctionBegin;
379:   if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS);
380:   PetscCall(PetscBLASIntCast(m, &bm));
381:   PetscCall(PetscBLASIntCast(n, &bn));
382: #if defined(PETSC_USE_COMPLEX)
383:   PetscCall(PetscMalloc2(m * n, &Js, m * n, &Jinvs));
384:   for (i = 0; i < m * n; i++) Js[i] = J[i];
385: #else
386:   Js    = (PetscReal *)J;
387:   Jinvs = Jinv;
388: #endif
389:   if (m == n) {
390:     PetscBLASInt *pivots;
391:     PetscScalar  *W;

393:     PetscCall(PetscMalloc2(m, &pivots, m, &W));

395:     PetscCall(PetscArraycpy(Jinvs, Js, m * m));
396:     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info));
397:     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscBLASInt_FMT, info);
398:     PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info));
399:     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscBLASInt_FMT, info);
400:     PetscCall(PetscFree2(pivots, W));
401:   } else if (m < n) {
402:     PetscScalar  *JJT;
403:     PetscBLASInt *pivots;
404:     PetscScalar  *W;

406:     PetscCall(PetscMalloc1(m * m, &JJT));
407:     PetscCall(PetscMalloc2(m, &pivots, m, &W));
408:     for (i = 0; i < m; i++) {
409:       for (j = 0; j < m; j++) {
410:         PetscScalar val = 0.;

412:         for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k];
413:         JJT[i * m + j] = val;
414:       }
415:     }

417:     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info));
418:     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscBLASInt_FMT, info);
419:     PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info));
420:     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscBLASInt_FMT, info);
421:     for (i = 0; i < n; i++) {
422:       for (j = 0; j < m; j++) {
423:         PetscScalar val = 0.;

425:         for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j];
426:         Jinvs[i * m + j] = val;
427:       }
428:     }
429:     PetscCall(PetscFree2(pivots, W));
430:     PetscCall(PetscFree(JJT));
431:   } else {
432:     PetscScalar  *JTJ;
433:     PetscBLASInt *pivots;
434:     PetscScalar  *W;

436:     PetscCall(PetscMalloc1(n * n, &JTJ));
437:     PetscCall(PetscMalloc2(n, &pivots, n, &W));
438:     for (i = 0; i < n; i++) {
439:       for (j = 0; j < n; j++) {
440:         PetscScalar val = 0.;

442:         for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j];
443:         JTJ[i * n + j] = val;
444:       }
445:     }

447:     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bn, pivots, &info));
448:     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscBLASInt_FMT, info);
449:     PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info));
450:     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscBLASInt_FMT, info);
451:     for (i = 0; i < n; i++) {
452:       for (j = 0; j < m; j++) {
453:         PetscScalar val = 0.;

455:         for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k];
456:         Jinvs[i * m + j] = val;
457:       }
458:     }
459:     PetscCall(PetscFree2(pivots, W));
460:     PetscCall(PetscFree(JTJ));
461:   }
462: #if defined(PETSC_USE_COMPLEX)
463:   for (i = 0; i < m * n; i++) Jinv[i] = PetscRealPart(Jinvs[i]);
464:   PetscCall(PetscFree2(Js, Jinvs));
465: #endif
466:   PetscFunctionReturn(PETSC_SUCCESS);
467: }

469: /*@
470:   PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation.

472:   Collective

474:   Input Parameters:
475: + q           - the quadrature functional
476: . imageDim    - the dimension of the image of the transformation
477: . origin      - a point in the original space
478: . originImage - the image of the origin under the transformation
479: . J           - the Jacobian of the image: an [imageDim x dim] matrix in row major order
480: - formDegree  - transform the quadrature weights as k-forms of this form degree (if the number of components is a multiple of (dim choose `formDegree`),
481:                 it is assumed that they represent multiple k-forms) [see `PetscDTAltVPullback()` for interpretation of `formDegree`]

483:   Output Parameter:
484: . Jinvstarq - a quadrature rule where each point is the image of a point in the original quadrature rule, and where the k-form weights have
485:   been pulled-back by the pseudoinverse of `J` to the k-form weights in the image space.

487:   Level: intermediate

489:   Note:
490:   The new quadrature rule will have a different number of components if spaces have different dimensions.
491:   For example, pushing a 2-form forward from a two dimensional space to a three dimensional space changes the number of components from 1 to 3.

493: .seealso: `PetscQuadrature`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
494: @*/
495: PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formDegree, PetscQuadrature *Jinvstarq)
496: {
497:   PetscInt         dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c;
498:   const PetscReal *points;
499:   const PetscReal *weights;
500:   PetscReal       *imagePoints, *imageWeights;
501:   PetscReal       *Jinv;
502:   PetscReal       *Jinvstar;

504:   PetscFunctionBegin;
506:   PetscCheck(imageDim >= PetscAbsInt(formDegree), PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %" PetscInt_FMT "-form in %" PetscInt_FMT " dimensions", PetscAbsInt(formDegree), imageDim);
507:   PetscCall(PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights));
508:   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &formSize));
509:   PetscCheck(Nc % formSize == 0, PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Number of components %" PetscInt_FMT " is not a multiple of formSize %" PetscInt_FMT, Nc, formSize);
510:   Ncopies = Nc / formSize;
511:   PetscCall(PetscDTBinomialInt(imageDim, PetscAbsInt(formDegree), &imageFormSize));
512:   imageNc = Ncopies * imageFormSize;
513:   PetscCall(PetscMalloc1(Npoints * imageDim, &imagePoints));
514:   PetscCall(PetscMalloc1(Npoints * imageNc, &imageWeights));
515:   PetscCall(PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar));
516:   PetscCall(PetscDTJacobianInverse_Internal(imageDim, dim, J, Jinv));
517:   PetscCall(PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formDegree, Jinvstar));
518:   for (pt = 0; pt < Npoints; pt++) {
519:     const PetscReal *point      = PetscSafePointerPlusOffset(points, pt * dim);
520:     PetscReal       *imagePoint = &imagePoints[pt * imageDim];

522:     for (i = 0; i < imageDim; i++) {
523:       PetscReal val = originImage[i];

525:       for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]);
526:       imagePoint[i] = val;
527:     }
528:     for (c = 0; c < Ncopies; c++) {
529:       const PetscReal *form      = &weights[pt * Nc + c * formSize];
530:       PetscReal       *imageForm = &imageWeights[pt * imageNc + c * imageFormSize];

532:       for (i = 0; i < imageFormSize; i++) {
533:         PetscReal val = 0.;

535:         for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j];
536:         imageForm[i] = val;
537:       }
538:     }
539:   }
540:   PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq));
541:   PetscCall(PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights));
542:   PetscCall(PetscFree2(Jinv, Jinvstar));
543:   PetscFunctionReturn(PETSC_SUCCESS);
544: }

546: /*@C
547:   PetscQuadratureSetData - Sets the data defining the quadrature

549:   Not Collective

551:   Input Parameters:
552: + q       - The `PetscQuadrature` object
553: . dim     - The spatial dimension
554: . Nc      - The number of components
555: . npoints - The number of quadrature points
556: . points  - The coordinates of each quadrature point
557: - weights - The weight of each quadrature point

559:   Level: intermediate

561:   Note:
562:   `q` owns the references to points and weights, so they must be allocated using `PetscMalloc()` and the user should not free them.

564: .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`, `PetscQuadratureGetData()`
565: @*/
566: PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) PeNSS
567: {
568:   PetscFunctionBegin;
570:   if (dim >= 0) q->dim = dim;
571:   if (Nc >= 0) q->Nc = Nc;
572:   if (npoints >= 0) q->numPoints = npoints;
573:   if (points) {
574:     PetscAssertPointer(points, 5);
575:     q->points = points;
576:   }
577:   if (weights) {
578:     PetscAssertPointer(weights, 6);
579:     q->weights = weights;
580:   }
581:   PetscFunctionReturn(PETSC_SUCCESS);
582: }

584: static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v)
585: {
586:   PetscInt          q, d, c;
587:   PetscViewerFormat format;

589:   PetscFunctionBegin;
590:   if (quad->Nc > 1)
591:     PetscCall(PetscViewerASCIIPrintf(v, "Quadrature on a %s of order %" PetscInt_FMT " on %" PetscInt_FMT " points (dim %" PetscInt_FMT ") with %" PetscInt_FMT " components\n", DMPolytopeTypes[quad->ct], quad->order, quad->numPoints, quad->dim, quad->Nc));
592:   else PetscCall(PetscViewerASCIIPrintf(v, "Quadrature on a %s of order %" PetscInt_FMT " on %" PetscInt_FMT " points (dim %" PetscInt_FMT ")\n", DMPolytopeTypes[quad->ct], quad->order, quad->numPoints, quad->dim));
593:   PetscCall(PetscViewerGetFormat(v, &format));
594:   if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS);
595:   for (q = 0; q < quad->numPoints; ++q) {
596:     PetscCall(PetscViewerASCIIPrintf(v, "p%" PetscInt_FMT " (", q));
597:     PetscCall(PetscViewerASCIIUseTabs(v, PETSC_FALSE));
598:     for (d = 0; d < quad->dim; ++d) {
599:       if (d) PetscCall(PetscViewerASCIIPrintf(v, ", "));
600:       PetscCall(PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q * quad->dim + d]));
601:     }
602:     PetscCall(PetscViewerASCIIPrintf(v, ") "));
603:     if (quad->Nc > 1) PetscCall(PetscViewerASCIIPrintf(v, "w%" PetscInt_FMT " (", q));
604:     for (c = 0; c < quad->Nc; ++c) {
605:       if (c) PetscCall(PetscViewerASCIIPrintf(v, ", "));
606:       PetscCall(PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q * quad->Nc + c]));
607:     }
608:     if (quad->Nc > 1) PetscCall(PetscViewerASCIIPrintf(v, ")"));
609:     PetscCall(PetscViewerASCIIPrintf(v, "\n"));
610:     PetscCall(PetscViewerASCIIUseTabs(v, PETSC_TRUE));
611:   }
612:   PetscFunctionReturn(PETSC_SUCCESS);
613: }

615: /*@
616:   PetscQuadratureView - View a `PetscQuadrature` object

618:   Collective

620:   Input Parameters:
621: + quad   - The `PetscQuadrature` object
622: - viewer - The `PetscViewer` object

624:   Level: beginner

626: .seealso: `PetscQuadrature`, `PetscViewer`, `PetscQuadratureCreate()`, `PetscQuadratureGetData()`
627: @*/
628: PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
629: {
630:   PetscBool isascii;

632:   PetscFunctionBegin;
635:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)quad), &viewer));
636:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
637:   PetscCall(PetscViewerASCIIPushTab(viewer));
638:   if (isascii) PetscCall(PetscQuadratureView_Ascii(quad, viewer));
639:   PetscCall(PetscViewerASCIIPopTab(viewer));
640:   PetscFunctionReturn(PETSC_SUCCESS);
641: }

643: /*@C
644:   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement

646:   Not Collective; No Fortran Support

648:   Input Parameters:
649: + q              - The original `PetscQuadrature`
650: . numSubelements - The number of subelements the original element is divided into
651: . v0             - An array of the initial points for each subelement
652: - jac            - An array of the Jacobian mappings from the reference to each subelement

654:   Output Parameter:
655: . qref - The dimension

657:   Level: intermediate

659:   Note:
660:   Together `v0` and `jac` define an affine mapping from the original reference element to each subelement

662: .seealso: `PetscQuadrature`, `PetscFECreate()`, `PetscSpaceGetDimension()`, `PetscDualSpaceGetDimension()`
663: @*/
664: PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
665: {
666:   DMPolytopeType   ct;
667:   const PetscReal *points, *weights;
668:   PetscReal       *pointsRef, *weightsRef;
669:   PetscInt         dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;

671:   PetscFunctionBegin;
673:   PetscAssertPointer(v0, 3);
674:   PetscAssertPointer(jac, 4);
675:   PetscAssertPointer(qref, 5);
676:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, qref));
677:   PetscCall(PetscQuadratureGetCellType(q, &ct));
678:   PetscCall(PetscQuadratureGetOrder(q, &order));
679:   PetscCall(PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights));
680:   npointsRef = npoints * numSubelements;
681:   PetscCall(PetscMalloc1(npointsRef * dim, &pointsRef));
682:   PetscCall(PetscMalloc1(npointsRef * Nc, &weightsRef));
683:   for (c = 0; c < numSubelements; ++c) {
684:     for (p = 0; p < npoints; ++p) {
685:       for (d = 0; d < dim; ++d) {
686:         pointsRef[(c * npoints + p) * dim + d] = v0[c * dim + d];
687:         for (e = 0; e < dim; ++e) pointsRef[(c * npoints + p) * dim + d] += jac[(c * dim + d) * dim + e] * (points[p * dim + e] + 1.0);
688:       }
689:       /* Could also use detJ here */
690:       for (cp = 0; cp < Nc; ++cp) weightsRef[(c * npoints + p) * Nc + cp] = weights[p * Nc + cp] / numSubelements;
691:     }
692:   }
693:   PetscCall(PetscQuadratureSetCellType(*qref, ct));
694:   PetscCall(PetscQuadratureSetOrder(*qref, order));
695:   PetscCall(PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef));
696:   PetscFunctionReturn(PETSC_SUCCESS);
697: }

699: /* Compute the coefficients for the Jacobi polynomial recurrence,

701:    J^{a,b}_n(x) = (cnm1 + cnm1x * x) * J^{a,b}_{n-1}(x) - cnm2 * J^{a,b}_{n-2}(x).
702:  */
703: #define PetscDTJacobiRecurrence_Internal(n, a, b, cnm1, cnm1x, cnm2) \
704:   do { \
705:     PetscReal _a = (a); \
706:     PetscReal _b = (b); \
707:     PetscReal _n = (n); \
708:     if (n == 1) { \
709:       (cnm1)  = (_a - _b) * 0.5; \
710:       (cnm1x) = (_a + _b + 2.) * 0.5; \
711:       (cnm2)  = 0.; \
712:     } else { \
713:       PetscReal _2n  = _n + _n; \
714:       PetscReal _d   = (_2n * (_n + _a + _b) * (_2n + _a + _b - 2)); \
715:       PetscReal _n1  = (_2n + _a + _b - 1.) * (_a * _a - _b * _b); \
716:       PetscReal _n1x = (_2n + _a + _b - 1.) * (_2n + _a + _b) * (_2n + _a + _b - 2); \
717:       PetscReal _n2  = 2. * ((_n + _a - 1.) * (_n + _b - 1.) * (_2n + _a + _b)); \
718:       (cnm1)         = _n1 / _d; \
719:       (cnm1x)        = _n1x / _d; \
720:       (cnm2)         = _n2 / _d; \
721:     } \
722:   } while (0)

724: /*@
725:   PetscDTJacobiNorm - Compute the weighted L2 norm of a Jacobi polynomial.

727:   $\| P^{\alpha,\beta}_n \|_{\alpha,\beta}^2 = \int_{-1}^1 (1 + x)^{\alpha} (1 - x)^{\beta} P^{\alpha,\beta}_n (x)^2 dx.$

729:   Input Parameters:
730: + alpha - the left exponent > -1
731: . beta  - the right exponent > -1
732: - n     - the polynomial degree

734:   Output Parameter:
735: . norm - the weighted L2 norm

737:   Level: beginner

739: .seealso: `PetscQuadrature`, `PetscDTJacobiEval()`
740: @*/
741: PetscErrorCode PetscDTJacobiNorm(PetscReal alpha, PetscReal beta, PetscInt n, PetscReal *norm)
742: {
743:   PetscReal twoab1;
744:   PetscReal gr;

746:   PetscFunctionBegin;
747:   PetscCheck(alpha > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exponent alpha %g <= -1. invalid", (double)alpha);
748:   PetscCheck(beta > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exponent beta %g <= -1. invalid", (double)beta);
749:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "n %" PetscInt_FMT " < 0 invalid", n);
750:   twoab1 = PetscPowReal(2., alpha + beta + 1.);
751: #if defined(PETSC_HAVE_LGAMMA)
752:   if (!n) {
753:     gr = PetscExpReal(PetscLGamma(alpha + 1.) + PetscLGamma(beta + 1.) - PetscLGamma(alpha + beta + 2.));
754:   } else {
755:     gr = PetscExpReal(PetscLGamma(n + alpha + 1.) + PetscLGamma(n + beta + 1.) - (PetscLGamma(n + 1.) + PetscLGamma(n + alpha + beta + 1.))) / (n + n + alpha + beta + 1.);
756:   }
757: #else
758:   {
759:     PetscInt alphai = (PetscInt)alpha;
760:     PetscInt betai  = (PetscInt)beta;
761:     PetscInt i;

763:     gr = n ? (1. / (n + n + alpha + beta + 1.)) : 1.;
764:     if ((PetscReal)alphai == alpha) {
765:       if (!n) {
766:         for (i = 0; i < alphai; i++) gr *= (i + 1.) / (beta + i + 1.);
767:         gr /= (alpha + beta + 1.);
768:       } else {
769:         for (i = 0; i < alphai; i++) gr *= (n + i + 1.) / (n + beta + i + 1.);
770:       }
771:     } else if ((PetscReal)betai == beta) {
772:       if (!n) {
773:         for (i = 0; i < betai; i++) gr *= (i + 1.) / (alpha + i + 2.);
774:         gr /= (alpha + beta + 1.);
775:       } else {
776:         for (i = 0; i < betai; i++) gr *= (n + i + 1.) / (n + alpha + i + 1.);
777:       }
778:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "lgamma() - math routine is unavailable.");
779:   }
780: #endif
781:   *norm = PetscSqrtReal(twoab1 * gr);
782:   PetscFunctionReturn(PETSC_SUCCESS);
783: }

785: static PetscErrorCode PetscDTJacobiEval_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscInt k, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *p)
786: {
787:   PetscReal ak, bk;
788:   PetscReal abk1;
789:   PetscInt  i, l, maxdegree;

791:   PetscFunctionBegin;
792:   maxdegree = degrees[ndegree - 1] - k;
793:   ak        = a + k;
794:   bk        = b + k;
795:   abk1      = a + b + k + 1.;
796:   if (maxdegree < 0) {
797:     for (i = 0; i < npoints; i++)
798:       for (l = 0; l < ndegree; l++) p[i * ndegree + l] = 0.;
799:     PetscFunctionReturn(PETSC_SUCCESS);
800:   }
801:   for (i = 0; i < npoints; i++) {
802:     PetscReal pm1, pm2, x;
803:     PetscReal cnm1, cnm1x, cnm2;
804:     PetscInt  j;

806:     x   = points[i];
807:     pm2 = 1.;
808:     PetscDTJacobiRecurrence_Internal(1, ak, bk, cnm1, cnm1x, cnm2);
809:     pm1 = (cnm1 + cnm1x * x);
810:     l   = 0;
811:     while (l < ndegree && degrees[l] - k < 0) p[l++] = 0.;
812:     while (l < ndegree && degrees[l] - k == 0) {
813:       p[l] = pm2;
814:       for (PetscInt m = 0; m < k; m++) p[l] *= (abk1 + m) * 0.5;
815:       l++;
816:     }
817:     while (l < ndegree && degrees[l] - k == 1) {
818:       p[l] = pm1;
819:       for (PetscInt m = 0; m < k; m++) p[l] *= (abk1 + 1 + m) * 0.5;
820:       l++;
821:     }
822:     for (j = 2; j <= maxdegree; j++) {
823:       PetscReal pp;

825:       PetscDTJacobiRecurrence_Internal(j, ak, bk, cnm1, cnm1x, cnm2);
826:       pp  = (cnm1 + cnm1x * x) * pm1 - cnm2 * pm2;
827:       pm2 = pm1;
828:       pm1 = pp;
829:       while (l < ndegree && degrees[l] - k == j) {
830:         p[l] = pp;
831:         for (PetscInt m = 0; m < k; m++) p[l] *= (abk1 + j + m) * 0.5;
832:         l++;
833:       }
834:     }
835:     p += ndegree;
836:   }
837:   PetscFunctionReturn(PETSC_SUCCESS);
838: }

840: /*@
841:   PetscDTJacobiEvalJet - Evaluate the jet (function and derivatives) of the Jacobi polynomials basis up to a given degree.

843:   Input Parameters:
844: + alpha   - the left exponent of the weight
845: . beta    - the right exponetn of the weight
846: . npoints - the number of points to evaluate the polynomials at
847: . points  - [npoints] array of point coordinates
848: . degree  - the maximm degree polynomial space to evaluate, (degree + 1) will be evaluated total.
849: - k       - the maximum derivative to evaluate in the jet, (k + 1) will be evaluated total.

851:   Output Parameter:
852: . p - an array containing the evaluations of the Jacobi polynomials's jets on the points.  the size is (degree + 1) x
853:       (k + 1) x npoints, which also describes the order of the dimensions of this three-dimensional array: the first
854:       (slowest varying) dimension is polynomial degree; the second dimension is derivative order; the third (fastest
855:       varying) dimension is the index of the evaluation point.

857:   Level: advanced

859:   Notes:
860:   The Jacobi polynomials with indices $\alpha$ and $\beta$ are orthogonal with respect to the
861:   weighted inner product $\langle f, g \rangle = \int_{-1}^1 (1+x)^{\alpha} (1-x)^{\beta} f(x)
862:   g(x) dx$.

864: .seealso: `PetscDTJacobiEval()`, `PetscDTPKDEvalJet()`
865: @*/
866: PetscErrorCode PetscDTJacobiEvalJet(PetscReal alpha, PetscReal beta, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt k, PetscReal p[])
867: {
868:   PetscInt   i, j, l;
869:   PetscInt  *degrees;
870:   PetscReal *psingle;

872:   PetscFunctionBegin;
873:   if (degree == 0) {
874:     PetscInt zero = 0;

876:     for (i = 0; i <= k; i++) PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, i, points, 1, &zero, &p[i * npoints]));
877:     PetscFunctionReturn(PETSC_SUCCESS);
878:   }
879:   PetscCall(PetscMalloc1(degree + 1, &degrees));
880:   PetscCall(PetscMalloc1((degree + 1) * npoints, &psingle));
881:   for (i = 0; i <= degree; i++) degrees[i] = i;
882:   for (i = 0; i <= k; i++) {
883:     PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, i, points, degree + 1, degrees, psingle));
884:     for (j = 0; j <= degree; j++) {
885:       for (l = 0; l < npoints; l++) p[(j * (k + 1) + i) * npoints + l] = psingle[l * (degree + 1) + j];
886:     }
887:   }
888:   PetscCall(PetscFree(psingle));
889:   PetscCall(PetscFree(degrees));
890:   PetscFunctionReturn(PETSC_SUCCESS);
891: }

893: /*@
894:   PetscDTJacobiEval - evaluate Jacobi polynomials for the weight function $(1.+x)^{\alpha} (1.-x)^{\beta}$ at a set of points
895:   at points

897:   Not Collective

899:   Input Parameters:
900: + npoints - number of spatial points to evaluate at
901: . alpha   - the left exponent > -1
902: . beta    - the right exponent > -1
903: . points  - array of locations to evaluate at
904: . ndegree - number of basis degrees to evaluate
905: - degrees - sorted array of degrees to evaluate

907:   Output Parameters:
908: + B  - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or `NULL`)
909: . D  - row-oriented derivative evaluation matrix (or `NULL`)
910: - D2 - row-oriented second derivative evaluation matrix (or `NULL`)

912:   Level: intermediate

914: .seealso: `PetscDTGaussQuadrature()`, `PetscDTLegendreEval()`
915: @*/
916: PetscErrorCode PetscDTJacobiEval(PetscInt npoints, PetscReal alpha, PetscReal beta, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PeOp PetscReal B[], PeOp PetscReal D[], PeOp PetscReal D2[])
917: {
918:   PetscFunctionBegin;
919:   PetscCheck(alpha > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "alpha must be > -1.");
920:   PetscCheck(beta > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "beta must be > -1.");
921:   if (!npoints || !ndegree) PetscFunctionReturn(PETSC_SUCCESS);
922:   if (B) PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, 0, points, ndegree, degrees, B));
923:   if (D) PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, 1, points, ndegree, degrees, D));
924:   if (D2) PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, 2, points, ndegree, degrees, D2));
925:   PetscFunctionReturn(PETSC_SUCCESS);
926: }

928: /*@
929:   PetscDTLegendreEval - evaluate Legendre polynomials at points

931:   Not Collective

933:   Input Parameters:
934: + npoints - number of spatial points to evaluate at
935: . points  - array of locations to evaluate at
936: . ndegree - number of basis degrees to evaluate
937: - degrees - sorted array of degrees to evaluate

939:   Output Parameters:
940: + B  - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension `npoints`*`ndegrees`, allocated by caller) (or `NULL`)
941: . D  - row-oriented derivative evaluation matrix (or `NULL`)
942: - D2 - row-oriented second derivative evaluation matrix (or `NULL`)

944:   Level: intermediate

946: .seealso: `PetscDTGaussQuadrature()`
947: @*/
948: PetscErrorCode PetscDTLegendreEval(PetscInt npoints, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PeOp PetscReal B[], PeOp PetscReal D[], PeOp PetscReal D2[])
949: {
950:   PetscFunctionBegin;
951:   PetscCall(PetscDTJacobiEval(npoints, 0., 0., points, ndegree, degrees, B, D, D2));
952:   PetscFunctionReturn(PETSC_SUCCESS);
953: }

955: /*@
956:   PetscDTIndexToGradedOrder - convert an index into a tuple of monomial degrees in a graded order (that is, if the degree sum of tuple x is less than the degree sum of tuple y,
957:   then the index of x is smaller than the index of y)

959:   Input Parameters:
960: + len   - the desired length of the degree tuple
961: - index - the index to convert: should be >= 0

963:   Output Parameter:
964: . degtup - filled with a tuple of degrees

966:   Level: beginner

968:   Note:
969:   For two tuples x and y with the same degree sum, partial degree sums over the final elements of the tuples
970:   acts as a tiebreaker.  For example, (2, 1, 1) and (1, 2, 1) have the same degree sum, but the degree sum over the
971:   last two elements is smaller for the former, so (2, 1, 1) < (1, 2, 1).

973: .seealso: `PetscDTGradedOrderToIndex()`
974: @*/
975: PetscErrorCode PetscDTIndexToGradedOrder(PetscInt len, PetscInt index, PetscInt degtup[])
976: {
977:   PetscInt i, total;
978:   PetscInt sum;

980:   PetscFunctionBeginHot;
981:   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
982:   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index must be non-negative");
983:   total = 1;
984:   sum   = 0;
985:   while (index >= total) {
986:     index -= total;
987:     total = (total * (len + sum)) / (sum + 1);
988:     sum++;
989:   }
990:   for (i = 0; i < len; i++) {
991:     PetscInt c;

993:     degtup[i] = sum;
994:     for (c = 0, total = 1; c < sum; c++) {
995:       /* going into the loop, total is the number of way to have a tuple of sum exactly c with length len - 1 - i */
996:       if (index < total) break;
997:       index -= total;
998:       total = (total * (len - 1 - i + c)) / (c + 1);
999:       degtup[i]--;
1000:     }
1001:     sum -= degtup[i];
1002:   }
1003:   PetscFunctionReturn(PETSC_SUCCESS);
1004: }

1006: /*@
1007:   PetscDTGradedOrderToIndex - convert a tuple into an index in a graded order, the inverse of `PetscDTIndexToGradedOrder()`.

1009:   Input Parameters:
1010: + len    - the length of the degree tuple
1011: - degtup - tuple with this length

1013:   Output Parameter:
1014: . index - index in graded order: >= 0

1016:   Level: beginner

1018:   Note:
1019:   For two tuples x and y with the same degree sum, partial degree sums over the final elements of the tuples
1020:   acts as a tiebreaker.  For example, (2, 1, 1) and (1, 2, 1) have the same degree sum, but the degree sum over the
1021:   last two elements is smaller for the former, so (2, 1, 1) < (1, 2, 1).

1023: .seealso: `PetscDTIndexToGradedOrder()`
1024: @*/
1025: PetscErrorCode PetscDTGradedOrderToIndex(PetscInt len, const PetscInt degtup[], PetscInt *index)
1026: {
1027:   PetscInt i, idx, sum, total;

1029:   PetscFunctionBeginHot;
1030:   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
1031:   for (i = 0, sum = 0; i < len; i++) sum += degtup[i];
1032:   idx   = 0;
1033:   total = 1;
1034:   for (i = 0; i < sum; i++) {
1035:     idx += total;
1036:     total = (total * (len + i)) / (i + 1);
1037:   }
1038:   for (i = 0; i < len - 1; i++) {
1039:     total = 1;
1040:     sum -= degtup[i];
1041:     for (PetscInt c = 0; c < sum; c++) {
1042:       idx += total;
1043:       total = (total * (len - 1 - i + c)) / (c + 1);
1044:     }
1045:   }
1046:   *index = idx;
1047:   PetscFunctionReturn(PETSC_SUCCESS);
1048: }

1050: static PetscBool PKDCite       = PETSC_FALSE;
1051: const char       PKDCitation[] = "@article{Kirby2010,\n"
1052:                                  "  title={Singularity-free evaluation of collapsed-coordinate orthogonal polynomials},\n"
1053:                                  "  author={Kirby, Robert C},\n"
1054:                                  "  journal={ACM Transactions on Mathematical Software (TOMS)},\n"
1055:                                  "  volume={37},\n"
1056:                                  "  number={1},\n"
1057:                                  "  pages={1--16},\n"
1058:                                  "  year={2010},\n"
1059:                                  "  publisher={ACM New York, NY, USA}\n}\n";

1061: /*@
1062:   PetscDTPKDEvalJet - Evaluate the jet (function and derivatives) of the Proriol-Koornwinder-Dubiner (PKD) basis for
1063:   the space of polynomials up to a given degree.

1065:   Input Parameters:
1066: + dim     - the number of variables in the multivariate polynomials
1067: . npoints - the number of points to evaluate the polynomials at
1068: . points  - [npoints x dim] array of point coordinates
1069: . degree  - the degree (sum of degrees on the variables in a monomial) of the polynomial space to evaluate.  There are ((dim + degree) choose dim) polynomials in this space.
1070: - k       - the maximum order partial derivative to evaluate in the jet.  There are (dim + k choose dim) partial derivatives
1071:             in the jet.  Choosing k = 0 means to evaluate just the function and no derivatives

1073:   Output Parameter:
1074: . p - an array containing the evaluations of the PKD polynomials' jets on the points.  The size is ((dim + degree)
1075:       choose dim) x ((dim + k) choose dim) x npoints, which also describes the order of the dimensions of this
1076:       three-dimensional array: the first (slowest varying) dimension is basis function index; the second dimension is jet
1077:       index; the third (fastest varying) dimension is the index of the evaluation point.

1079:   Level: advanced

1081:   Notes:
1082:   The PKD basis is L2-orthonormal on the biunit simplex (which is used as the reference element
1083:   for finite elements in PETSc), which makes it a stable basis to use for evaluating
1084:   polynomials in that domain.

1086:   The ordering of the basis functions, and the ordering of the derivatives in the jet, both follow the graded
1087:   ordering of `PetscDTIndexToGradedOrder()` and `PetscDTGradedOrderToIndex()`.  For example, in 3D, the polynomial with
1088:   leading monomial x^2,y^0,z^1, which has degree tuple (2,0,1), which by `PetscDTGradedOrderToIndex()` has index 12 (it is the 13th basis function in the space);
1089:   the partial derivative $\partial_x \partial_z$ has order tuple (1,0,1), appears at index 6 in the jet (it is the 7th partial derivative in the jet).

1091:   The implementation uses Kirby's singularity-free evaluation algorithm, <https://doi.org/10.1145/1644001.1644006>.

1093: .seealso: `PetscDTGradedOrderToIndex()`, `PetscDTIndexToGradedOrder()`, `PetscDTJacobiEvalJet()`
1094: @*/
1095: PetscErrorCode PetscDTPKDEvalJet(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt k, PetscReal p[])
1096: {
1097:   PetscInt   degidx, kidx, d, pt;
1098:   PetscInt   Nk, Ndeg;
1099:   PetscInt  *ktup, *degtup;
1100:   PetscReal *scales, initscale, scaleexp;

1102:   PetscFunctionBegin;
1103:   PetscCall(PetscCitationsRegister(PKDCitation, &PKDCite));
1104:   PetscCall(PetscDTBinomialInt(dim + k, k, &Nk));
1105:   PetscCall(PetscDTBinomialInt(degree + dim, degree, &Ndeg));
1106:   PetscCall(PetscMalloc2(dim, &degtup, dim, &ktup));
1107:   PetscCall(PetscMalloc1(Ndeg, &scales));
1108:   initscale = 1.;
1109:   if (dim > 1) {
1110:     PetscCall(PetscDTBinomial(dim, 2, &scaleexp));
1111:     initscale = PetscPowReal(2., scaleexp * 0.5);
1112:   }
1113:   for (degidx = 0; degidx < Ndeg; degidx++) {
1114:     PetscInt  e;
1115:     PetscInt  m1idx = -1, m2idx = -1;
1116:     PetscInt  n;
1117:     PetscInt  degsum;
1118:     PetscReal alpha;
1119:     PetscReal cnm1, cnm1x, cnm2;
1120:     PetscReal norm;

1122:     PetscCall(PetscDTIndexToGradedOrder(dim, degidx, degtup));
1123:     for (d = dim - 1; d >= 0; d--)
1124:       if (degtup[d]) break;
1125:     if (d < 0) { /* constant is 1 everywhere, all derivatives are zero */
1126:       scales[degidx] = initscale;
1127:       for (e = 0; e < dim; e++) {
1128:         PetscCall(PetscDTJacobiNorm(e, 0., 0, &norm));
1129:         scales[degidx] /= norm;
1130:       }
1131:       for (PetscInt i = 0; i < npoints; i++) p[degidx * Nk * npoints + i] = 1.;
1132:       for (PetscInt i = 0; i < (Nk - 1) * npoints; i++) p[(degidx * Nk + 1) * npoints + i] = 0.;
1133:       continue;
1134:     }
1135:     n = degtup[d];
1136:     degtup[d]--;
1137:     PetscCall(PetscDTGradedOrderToIndex(dim, degtup, &m1idx));
1138:     if (degtup[d] > 0) {
1139:       degtup[d]--;
1140:       PetscCall(PetscDTGradedOrderToIndex(dim, degtup, &m2idx));
1141:       degtup[d]++;
1142:     }
1143:     degtup[d]++;
1144:     for (e = 0, degsum = 0; e < d; e++) degsum += degtup[e];
1145:     alpha = 2 * degsum + d;
1146:     PetscDTJacobiRecurrence_Internal(n, alpha, 0., cnm1, cnm1x, cnm2);

1148:     scales[degidx] = initscale;
1149:     for (e = 0, degsum = 0; e < dim; e++) {
1150:       PetscReal ealpha;
1151:       PetscReal enorm;

1153:       ealpha = 2 * degsum + e;
1154:       for (PetscInt f = 0; f < degsum; f++) scales[degidx] *= 2.;
1155:       PetscCall(PetscDTJacobiNorm(ealpha, 0., degtup[e], &enorm));
1156:       scales[degidx] /= enorm;
1157:       degsum += degtup[e];
1158:     }

1160:     for (pt = 0; pt < npoints; pt++) {
1161:       /* compute the multipliers */
1162:       PetscReal thetanm1, thetanm1x, thetanm2;

1164:       thetanm1x = dim - (d + 1) + 2. * points[pt * dim + d];
1165:       for (e = d + 1; e < dim; e++) thetanm1x += points[pt * dim + e];
1166:       thetanm1x *= 0.5;
1167:       thetanm1 = (2. - (dim - (d + 1)));
1168:       for (e = d + 1; e < dim; e++) thetanm1 -= points[pt * dim + e];
1169:       thetanm1 *= 0.5;
1170:       thetanm2 = thetanm1 * thetanm1;

1172:       for (kidx = 0; kidx < Nk; kidx++) {
1173:         PetscCall(PetscDTIndexToGradedOrder(dim, kidx, ktup));
1174:         /* first sum in the same derivative terms */
1175:         p[(degidx * Nk + kidx) * npoints + pt] = (cnm1 * thetanm1 + cnm1x * thetanm1x) * p[(m1idx * Nk + kidx) * npoints + pt];
1176:         if (m2idx >= 0) p[(degidx * Nk + kidx) * npoints + pt] -= cnm2 * thetanm2 * p[(m2idx * Nk + kidx) * npoints + pt];

1178:         for (PetscInt f = d; f < dim; f++) {
1179:           PetscInt km1idx, mplty = ktup[f];

1181:           if (!mplty) continue;
1182:           ktup[f]--;
1183:           PetscCall(PetscDTGradedOrderToIndex(dim, ktup, &km1idx));

1185:           /* the derivative of  cnm1x * thetanm1x  wrt x variable f is 0.5 * cnm1x if f > d otherwise it is cnm1x */
1186:           /* the derivative of  cnm1  * thetanm1   wrt x variable f is 0 if f == d, otherwise it is -0.5 * cnm1 */
1187:           /* the derivative of -cnm2  * thetanm2   wrt x variable f is 0 if f == d, otherwise it is cnm2 * thetanm1 */
1188:           if (f > d) {
1189:             p[(degidx * Nk + kidx) * npoints + pt] += mplty * 0.5 * (cnm1x - cnm1) * p[(m1idx * Nk + km1idx) * npoints + pt];
1190:             if (m2idx >= 0) {
1191:               p[(degidx * Nk + kidx) * npoints + pt] += mplty * cnm2 * thetanm1 * p[(m2idx * Nk + km1idx) * npoints + pt];
1192:               /* second derivatives of -cnm2  * thetanm2   wrt x variable f,f2 is like - 0.5 * cnm2 */
1193:               for (PetscInt f2 = f; f2 < dim; f2++) {
1194:                 PetscInt km2idx, mplty2 = ktup[f2];
1195:                 PetscInt factor;

1197:                 if (!mplty2) continue;
1198:                 ktup[f2]--;
1199:                 PetscCall(PetscDTGradedOrderToIndex(dim, ktup, &km2idx));

1201:                 factor = mplty * mplty2;
1202:                 if (f == f2) factor /= 2;
1203:                 p[(degidx * Nk + kidx) * npoints + pt] -= 0.5 * factor * cnm2 * p[(m2idx * Nk + km2idx) * npoints + pt];
1204:                 ktup[f2]++;
1205:               }
1206:             }
1207:           } else {
1208:             p[(degidx * Nk + kidx) * npoints + pt] += mplty * cnm1x * p[(m1idx * Nk + km1idx) * npoints + pt];
1209:           }
1210:           ktup[f]++;
1211:         }
1212:       }
1213:     }
1214:   }
1215:   for (degidx = 0; degidx < Ndeg; degidx++) {
1216:     PetscReal scale = scales[degidx];

1218:     for (PetscInt i = 0; i < Nk * npoints; i++) p[degidx * Nk * npoints + i] *= scale;
1219:   }
1220:   PetscCall(PetscFree(scales));
1221:   PetscCall(PetscFree2(degtup, ktup));
1222:   PetscFunctionReturn(PETSC_SUCCESS);
1223: }

1225: /*@
1226:   PetscDTPTrimmedSize - The size of the trimmed polynomial space of k-forms with a given degree and form degree,
1227:   which can be evaluated in `PetscDTPTrimmedEvalJet()`.

1229:   Input Parameters:
1230: + dim        - the number of variables in the multivariate polynomials
1231: . degree     - the degree (sum of degrees on the variables in a monomial) of the trimmed polynomial space.
1232: - formDegree - the degree of the form

1234:   Output Parameter:
1235: . size - The number ((`dim` + `degree`) choose (`dim` + `formDegree`)) x ((`degree` + `formDegree` - 1) choose (`formDegree`))

1237:   Level: advanced

1239: .seealso: `PetscDTPTrimmedEvalJet()`
1240: @*/
1241: PetscErrorCode PetscDTPTrimmedSize(PetscInt dim, PetscInt degree, PetscInt formDegree, PetscInt *size)
1242: {
1243:   PetscInt Nrk, Nbpt; // number of trimmed polynomials

1245:   PetscFunctionBegin;
1246:   formDegree = PetscAbsInt(formDegree);
1247:   PetscCall(PetscDTBinomialInt(degree + dim, degree + formDegree, &Nbpt));
1248:   PetscCall(PetscDTBinomialInt(degree + formDegree - 1, formDegree, &Nrk));
1249:   Nbpt *= Nrk;
1250:   *size = Nbpt;
1251:   PetscFunctionReturn(PETSC_SUCCESS);
1252: }

1254: /* there was a reference implementation based on section 4.4 of Arnold, Falk & Winther (acta numerica, 2006), but it
1255:  * was inferior to this implementation */
1256: static PetscErrorCode PetscDTPTrimmedEvalJet_Internal(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt formDegree, PetscInt jetDegree, PetscReal p[])
1257: {
1258:   PetscInt  formDegreeOrig = formDegree;
1259:   PetscBool formNegative   = (formDegreeOrig < 0) ? PETSC_TRUE : PETSC_FALSE;

1261:   PetscFunctionBegin;
1262:   formDegree = PetscAbsInt(formDegreeOrig);
1263:   if (formDegree == 0) {
1264:     PetscCall(PetscDTPKDEvalJet(dim, npoints, points, degree, jetDegree, p));
1265:     PetscFunctionReturn(PETSC_SUCCESS);
1266:   }
1267:   if (formDegree == dim) {
1268:     PetscCall(PetscDTPKDEvalJet(dim, npoints, points, degree - 1, jetDegree, p));
1269:     PetscFunctionReturn(PETSC_SUCCESS);
1270:   }
1271:   PetscInt Nbpt;
1272:   PetscCall(PetscDTPTrimmedSize(dim, degree, formDegree, &Nbpt));
1273:   PetscInt Nf;
1274:   PetscCall(PetscDTBinomialInt(dim, formDegree, &Nf));
1275:   PetscInt Nk;
1276:   PetscCall(PetscDTBinomialInt(dim + jetDegree, dim, &Nk));
1277:   PetscCall(PetscArrayzero(p, Nbpt * Nf * Nk * npoints));

1279:   PetscInt Nbpm1; // number of scalar polynomials up to degree - 1;
1280:   PetscCall(PetscDTBinomialInt(dim + degree - 1, dim, &Nbpm1));
1281:   PetscReal *p_scalar;
1282:   PetscCall(PetscMalloc1(Nbpm1 * Nk * npoints, &p_scalar));
1283:   PetscCall(PetscDTPKDEvalJet(dim, npoints, points, degree - 1, jetDegree, p_scalar));
1284:   PetscInt total = 0;
1285:   // First add the full polynomials up to degree - 1 into the basis: take the scalar
1286:   // and copy one for each form component
1287:   for (PetscInt i = 0; i < Nbpm1; i++) {
1288:     const PetscReal *src = &p_scalar[i * Nk * npoints];
1289:     for (PetscInt f = 0; f < Nf; f++) {
1290:       PetscReal *dest = &p[(total++ * Nf + f) * Nk * npoints];
1291:       PetscCall(PetscArraycpy(dest, src, Nk * npoints));
1292:     }
1293:   }
1294:   PetscInt *form_atoms;
1295:   PetscCall(PetscMalloc1(formDegree + 1, &form_atoms));
1296:   // construct the interior product pattern
1297:   PetscInt (*pattern)[3];
1298:   PetscInt Nf1; // number of formDegree + 1 forms
1299:   PetscCall(PetscDTBinomialInt(dim, formDegree + 1, &Nf1));
1300:   PetscInt nnz = Nf1 * (formDegree + 1);
1301:   PetscCall(PetscMalloc1(Nf1 * (formDegree + 1), &pattern));
1302:   PetscCall(PetscDTAltVInteriorPattern(dim, formDegree + 1, pattern));
1303:   PetscReal centroid = (1. - dim) / (dim + 1.);
1304:   PetscInt *deriv;
1305:   PetscCall(PetscMalloc1(dim, &deriv));
1306:   for (PetscInt d = dim; d >= formDegree + 1; d--) {
1307:     PetscInt Nfd1; // number of formDegree + 1 forms in dimension d that include dx_0
1308:                    // (equal to the number of formDegree forms in dimension d-1)
1309:     PetscCall(PetscDTBinomialInt(d - 1, formDegree, &Nfd1));
1310:     // The number of homogeneous (degree-1) scalar polynomials in d variables
1311:     PetscInt Nh;
1312:     PetscCall(PetscDTBinomialInt(d - 1 + degree - 1, d - 1, &Nh));
1313:     const PetscReal *h_scalar = &p_scalar[(Nbpm1 - Nh) * Nk * npoints];
1314:     for (PetscInt b = 0; b < Nh; b++) {
1315:       const PetscReal *h_s = &h_scalar[b * Nk * npoints];
1316:       for (PetscInt f = 0; f < Nfd1; f++) {
1317:         // construct all formDegree+1 forms that start with dx_(dim - d) /\ ...
1318:         form_atoms[0] = dim - d;
1319:         PetscCall(PetscDTEnumSubset(d - 1, formDegree, f, &form_atoms[1]));
1320:         for (PetscInt i = 0; i < formDegree; i++) form_atoms[1 + i] += form_atoms[0] + 1;
1321:         PetscInt f_ind; // index of the resulting form
1322:         PetscCall(PetscDTSubsetIndex(dim, formDegree + 1, form_atoms, &f_ind));
1323:         PetscReal *p_f = &p[total++ * Nf * Nk * npoints];
1324:         for (PetscInt nz = 0; nz < nnz; nz++) {
1325:           PetscInt  i     = pattern[nz][0]; // formDegree component
1326:           PetscInt  j     = pattern[nz][1]; // (formDegree + 1) component
1327:           PetscInt  v     = pattern[nz][2]; // coordinate component
1328:           PetscReal scale = v < 0 ? -1. : 1.;

1330:           i     = formNegative ? (Nf - 1 - i) : i;
1331:           scale = (formNegative && (i & 1)) ? -scale : scale;
1332:           v     = v < 0 ? -(v + 1) : v;
1333:           if (j != f_ind) continue;
1334:           PetscReal *p_i = &p_f[i * Nk * npoints];
1335:           for (PetscInt jet = 0; jet < Nk; jet++) {
1336:             const PetscReal *h_jet = &h_s[jet * npoints];
1337:             PetscReal       *p_jet = &p_i[jet * npoints];

1339:             for (PetscInt pt = 0; pt < npoints; pt++) p_jet[pt] += scale * h_jet[pt] * (points[pt * dim + v] - centroid);
1340:             PetscCall(PetscDTIndexToGradedOrder(dim, jet, deriv));
1341:             deriv[v]++;
1342:             PetscReal mult = deriv[v];
1343:             PetscInt  l;
1344:             PetscCall(PetscDTGradedOrderToIndex(dim, deriv, &l));
1345:             if (l >= Nk) continue;
1346:             p_jet = &p_i[l * npoints];
1347:             for (PetscInt pt = 0; pt < npoints; pt++) p_jet[pt] += scale * mult * h_jet[pt];
1348:             deriv[v]--;
1349:           }
1350:         }
1351:       }
1352:     }
1353:   }
1354:   PetscCheck(total == Nbpt, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrectly counted P trimmed polynomials");
1355:   PetscCall(PetscFree(deriv));
1356:   PetscCall(PetscFree(pattern));
1357:   PetscCall(PetscFree(form_atoms));
1358:   PetscCall(PetscFree(p_scalar));
1359:   PetscFunctionReturn(PETSC_SUCCESS);
1360: }

1362: /*@
1363:   PetscDTPTrimmedEvalJet - Evaluate the jet (function and derivatives) of a basis of the trimmed polynomial k-forms up to
1364:   a given degree.

1366:   Input Parameters:
1367: + dim        - the number of variables in the multivariate polynomials
1368: . npoints    - the number of points to evaluate the polynomials at
1369: . points     - [npoints x dim] array of point coordinates
1370: . degree     - the degree (sum of degrees on the variables in a monomial) of the trimmed polynomial space to evaluate.
1371:            There are ((dim + degree) choose (dim + formDegree)) x ((degree + formDegree - 1) choose (formDegree)) polynomials in this space.
1372:            (You can use `PetscDTPTrimmedSize()` to compute this size.)
1373: . formDegree - the degree of the form
1374: - jetDegree  - the maximum order partial derivative to evaluate in the jet.  There are ((dim + jetDegree) choose dim) partial derivatives
1375:               in the jet.  Choosing jetDegree = 0 means to evaluate just the function and no derivatives

1377:   Output Parameter:
1378: . p - an array containing the evaluations of the PKD polynomials' jets on the points.

1380:   Level: advanced

1382:   Notes:
1383:   The size of `p` is `PetscDTPTrimmedSize()` x ((dim + formDegree) choose dim) x ((dim + k)
1384:   choose dim) x npoints,which also describes the order of the dimensions of this
1385:   four-dimensional array\:

1387:   the first (slowest varying) dimension is basis function index;
1388:   the second dimension is component of the form;
1389:   the third dimension is jet index;
1390:   the fourth (fastest varying) dimension is the index of the evaluation point.

1392:   The ordering of the basis functions is not graded, so the basis functions are not nested by degree like `PetscDTPKDEvalJet()`.
1393:   The basis functions are not an L2-orthonormal basis on any particular domain.

1395:   The implementation is based on the description of the trimmed polynomials up to degree r as
1396:   the direct sum of polynomials up to degree (r-1) and the Koszul differential applied to
1397:   homogeneous polynomials of degree (r-1).

1399: .seealso: `PetscDTPKDEvalJet()`, `PetscDTPTrimmedSize()`
1400: @*/
1401: PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt formDegree, PetscInt jetDegree, PetscReal p[])
1402: {
1403:   PetscFunctionBegin;
1404:   PetscCall(PetscDTPTrimmedEvalJet_Internal(dim, npoints, points, degree, formDegree, jetDegree, p));
1405:   PetscFunctionReturn(PETSC_SUCCESS);
1406: }

1408: /* solve the symmetric tridiagonal eigenvalue system, writing the eigenvalues into eigs and the eigenvectors into V
1409:  * with lds n; diag and subdiag are overwritten */
1410: static PetscErrorCode PetscDTSymmetricTridiagonalEigensolve(PetscInt n, PetscReal diag[], PetscReal subdiag[], PetscReal eigs[], PetscScalar V[])
1411: {
1412:   char          jobz   = 'V'; /* eigenvalues and eigenvectors */
1413:   char          range  = 'A'; /* all eigenvalues will be found */
1414:   PetscReal     VL     = 0.;  /* ignored because range is 'A' */
1415:   PetscReal     VU     = 0.;  /* ignored because range is 'A' */
1416:   PetscBLASInt  IL     = 0;   /* ignored because range is 'A' */
1417:   PetscBLASInt  IU     = 0;   /* ignored because range is 'A' */
1418:   PetscReal     abstol = 0.;  /* unused */
1419:   PetscBLASInt  bn, bm, ldz;  /* bm will equal bn on exit */
1420:   PetscBLASInt *isuppz;
1421:   PetscBLASInt  lwork, liwork;
1422:   PetscReal     workquery;
1423:   PetscBLASInt  iworkquery;
1424:   PetscBLASInt *iwork;
1425:   PetscBLASInt  info;
1426:   PetscReal    *work = NULL;

1428:   PetscFunctionBegin;
1429: #if !defined(PETSCDTGAUSSIANQUADRATURE_EIG)
1430:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
1431: #endif
1432:   PetscCall(PetscBLASIntCast(n, &bn));
1433:   PetscCall(PetscBLASIntCast(n, &ldz));
1434: #if !defined(PETSC_MISSING_LAPACK_STEGR)
1435:   PetscCall(PetscMalloc1(2 * n, &isuppz));
1436:   lwork  = -1;
1437:   liwork = -1;
1438:   PetscCallBLAS("LAPACKstegr", LAPACKstegr_(&jobz, &range, &bn, diag, subdiag, &VL, &VU, &IL, &IU, &abstol, &bm, eigs, V, &ldz, isuppz, &workquery, &lwork, &iworkquery, &liwork, &info));
1439:   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xSTEGR error");
1440:   lwork  = (PetscBLASInt)workquery;
1441:   liwork = iworkquery;
1442:   PetscCall(PetscMalloc2(lwork, &work, liwork, &iwork));
1443:   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1444:   PetscCallBLAS("LAPACKstegr", LAPACKstegr_(&jobz, &range, &bn, diag, subdiag, &VL, &VU, &IL, &IU, &abstol, &bm, eigs, V, &ldz, isuppz, work, &lwork, iwork, &liwork, &info));
1445:   PetscCall(PetscFPTrapPop());
1446:   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xSTEGR error");
1447:   PetscCall(PetscFree2(work, iwork));
1448:   PetscCall(PetscFree(isuppz));
1449: #elif !defined(PETSC_MISSING_LAPACK_STEQR)
1450:   jobz = 'I'; /* Compute eigenvalues and eigenvectors of the
1451:                  tridiagonal matrix.  Z is initialized to the identity
1452:                  matrix. */
1453:   PetscCall(PetscMalloc1(PetscMax(1, 2 * n - 2), &work));
1454:   PetscCallBLAS("LAPACKsteqr", LAPACKsteqr_("I", &bn, diag, subdiag, V, &ldz, work, &info));
1455:   PetscCall(PetscFPTrapPop());
1456:   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xSTEQR error");
1457:   PetscCall(PetscFree(work));
1458:   PetscCall(PetscArraycpy(eigs, diag, n));
1459: #endif
1460:   PetscFunctionReturn(PETSC_SUCCESS);
1461: }

1463: /* Formula for the weights at the endpoints (-1 and 1) of Gauss-Lobatto-Jacobi
1464:  * quadrature rules on the interval [-1, 1] */
1465: static PetscErrorCode PetscDTGaussLobattoJacobiEndweights_Internal(PetscInt n, PetscReal alpha, PetscReal beta, PetscReal *leftw, PetscReal *rightw)
1466: {
1467:   PetscReal twoab1;
1468:   PetscInt  m = n - 2;
1469:   PetscReal a = alpha + 1.;
1470:   PetscReal b = beta + 1.;
1471:   PetscReal gra, grb;

1473:   PetscFunctionBegin;
1474:   twoab1 = PetscPowReal(2., a + b - 1.);
1475: #if defined(PETSC_HAVE_LGAMMA)
1476:   grb = PetscExpReal(2. * PetscLGamma(b + 1.) + PetscLGamma(m + 1.) + PetscLGamma(m + a + 1.) - (PetscLGamma(m + b + 1) + PetscLGamma(m + a + b + 1.)));
1477:   gra = PetscExpReal(2. * PetscLGamma(a + 1.) + PetscLGamma(m + 1.) + PetscLGamma(m + b + 1.) - (PetscLGamma(m + a + 1) + PetscLGamma(m + a + b + 1.)));
1478: #else
1479:   {
1480:     PetscReal binom1, binom2;
1481:     PetscInt  alphai = (PetscInt)alpha;
1482:     PetscInt  betai  = (PetscInt)beta;

1484:     PetscCheck((PetscReal)alphai == alpha && (PetscReal)betai == beta, PETSC_COMM_SELF, PETSC_ERR_SUP, "lgamma() - math routine is unavailable.");
1485:     PetscCall(PetscDTBinomial(m + b, b, &binom1));
1486:     PetscCall(PetscDTBinomial(m + a + b, b, &binom2));
1487:     grb = 1. / (binom1 * binom2);
1488:     PetscCall(PetscDTBinomial(m + a, a, &binom1));
1489:     PetscCall(PetscDTBinomial(m + a + b, a, &binom2));
1490:     gra = 1. / (binom1 * binom2);
1491:   }
1492: #endif
1493:   *leftw  = twoab1 * grb / b;
1494:   *rightw = twoab1 * gra / a;
1495:   PetscFunctionReturn(PETSC_SUCCESS);
1496: }

1498: /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
1499:    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
1500: static inline PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
1501: {
1502:   PetscReal pn1, pn2;
1503:   PetscReal cnm1, cnm1x, cnm2;

1505:   PetscFunctionBegin;
1506:   if (!n) {
1507:     *P = 1.0;
1508:     PetscFunctionReturn(PETSC_SUCCESS);
1509:   }
1510:   PetscDTJacobiRecurrence_Internal(1, a, b, cnm1, cnm1x, cnm2);
1511:   pn2 = 1.;
1512:   pn1 = cnm1 + cnm1x * x;
1513:   if (n == 1) {
1514:     *P = pn1;
1515:     PetscFunctionReturn(PETSC_SUCCESS);
1516:   }
1517:   *P = 0.0;
1518:   for (PetscInt k = 2; k < n + 1; ++k) {
1519:     PetscDTJacobiRecurrence_Internal(k, a, b, cnm1, cnm1x, cnm2);

1521:     *P  = (cnm1 + cnm1x * x) * pn1 - cnm2 * pn2;
1522:     pn2 = pn1;
1523:     pn1 = *P;
1524:   }
1525:   PetscFunctionReturn(PETSC_SUCCESS);
1526: }

1528: /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
1529: static inline PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscInt k, PetscReal *P)
1530: {
1531:   PetscReal nP;
1532:   PetscInt  i;

1534:   PetscFunctionBegin;
1535:   *P = 0.0;
1536:   if (k > n) PetscFunctionReturn(PETSC_SUCCESS);
1537:   PetscCall(PetscDTComputeJacobi(a + k, b + k, n - k, x, &nP));
1538:   for (i = 0; i < k; i++) nP *= (a + b + n + 1. + i) * 0.5;
1539:   *P = nP;
1540:   PetscFunctionReturn(PETSC_SUCCESS);
1541: }

1543: static PetscErrorCode PetscDTGaussJacobiQuadrature_Newton_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[])
1544: {
1545:   PetscInt  maxIter = 100;
1546:   PetscReal eps     = PetscExpReal(0.75 * PetscLogReal(PETSC_MACHINE_EPSILON));
1547:   PetscReal a1, a6, gf;

1549:   PetscFunctionBegin;
1550:   a1 = PetscPowReal(2.0, a + b + 1);
1551: #if defined(PETSC_HAVE_LGAMMA)
1552:   {
1553:     PetscReal a2, a3, a4, a5;
1554:     a2 = PetscLGamma(a + npoints + 1);
1555:     a3 = PetscLGamma(b + npoints + 1);
1556:     a4 = PetscLGamma(a + b + npoints + 1);
1557:     a5 = PetscLGamma(npoints + 1);
1558:     gf = PetscExpReal(a2 + a3 - (a4 + a5));
1559:   }
1560: #else
1561:   {
1562:     PetscInt ia, ib;

1564:     ia = (PetscInt)a;
1565:     ib = (PetscInt)b;
1566:     gf = 1.;
1567:     if (ia == a && ia >= 0) { /* compute ratio of rising factorals wrt a */
1568:       for (PetscInt k = 0; k < ia; k++) gf *= (npoints + 1. + k) / (npoints + b + 1. + k);
1569:     } else if (b == b && ib >= 0) { /* compute ratio of rising factorials wrt b */
1570:       for (PetscInt k = 0; k < ib; k++) gf *= (npoints + 1. + k) / (npoints + a + 1. + k);
1571:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "lgamma() - math routine is unavailable.");
1572:   }
1573: #endif

1575:   a6 = a1 * gf;
1576:   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
1577:    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
1578:   for (PetscInt k = 0; k < npoints; ++k) {
1579:     PetscReal r = PetscCosReal(PETSC_PI * (1. - (4. * k + 3. + 2. * b) / (4. * npoints + 2. * (a + b + 1.)))), dP;

1581:     if (k > 0) r = 0.5 * (r + x[k - 1]);
1582:     for (PetscInt j = 0; j < maxIter; ++j) {
1583:       PetscReal s = 0.0, delta, f, fp;
1584:       PetscInt  i;

1586:       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
1587:       PetscCall(PetscDTComputeJacobi(a, b, npoints, r, &f));
1588:       PetscCall(PetscDTComputeJacobiDerivative(a, b, npoints, r, 1, &fp));
1589:       delta = f / (fp - f * s);
1590:       r     = r - delta;
1591:       if (PetscAbsReal(delta) < eps) break;
1592:     }
1593:     x[k] = r;
1594:     PetscCall(PetscDTComputeJacobiDerivative(a, b, npoints, x[k], 1, &dP));
1595:     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
1596:   }
1597:   PetscFunctionReturn(PETSC_SUCCESS);
1598: }

1600: /* Compute the diagonals of the Jacobi matrix used in Golub & Welsch algorithms for Gauss-Jacobi
1601:  * quadrature weight calculations on [-1,1] for exponents (1. + x)^a (1.-x)^b */
1602: static PetscErrorCode PetscDTJacobiMatrix_Internal(PetscInt nPoints, PetscReal a, PetscReal b, PetscReal *d, PetscReal *s)
1603: {
1604:   PetscInt i;

1606:   PetscFunctionBegin;
1607:   for (i = 0; i < nPoints; i++) {
1608:     PetscReal A, B, C;

1610:     PetscDTJacobiRecurrence_Internal(i + 1, a, b, A, B, C);
1611:     d[i] = -A / B;
1612:     if (i) s[i - 1] *= C / B;
1613:     if (i < nPoints - 1) s[i] = 1. / B;
1614:   }
1615:   PetscFunctionReturn(PETSC_SUCCESS);
1616: }

1618: static PetscErrorCode PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
1619: {
1620:   PetscReal mu0;
1621:   PetscReal ga, gb, gab;
1622:   PetscInt  i;

1624:   PetscFunctionBegin;
1625:   PetscCall(PetscCitationsRegister(GolubWelschCitation, &GolubWelschCite));

1627: #if defined(PETSC_HAVE_TGAMMA)
1628:   ga  = PetscTGamma(a + 1);
1629:   gb  = PetscTGamma(b + 1);
1630:   gab = PetscTGamma(a + b + 2);
1631: #else
1632:   {
1633:     PetscInt ia = (PetscInt)a, ib = (PetscInt)b;

1635:     PetscCheck(ia == a && ib == b && ia + 1 > 0 && ib + 1 > 0 && ia + ib + 2 > 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "tgamma() - math routine is unavailable.");
1636:     /* All gamma(x) terms are (x-1)! terms */
1637:     PetscCall(PetscDTFactorial(ia, &ga));
1638:     PetscCall(PetscDTFactorial(ib, &gb));
1639:     PetscCall(PetscDTFactorial(ia + ib + 1, &gab));
1640:   }
1641: #endif
1642:   mu0 = PetscPowReal(2., a + b + 1.) * ga * gb / gab;

1644: #if defined(PETSCDTGAUSSIANQUADRATURE_EIG)
1645:   {
1646:     PetscReal   *diag, *subdiag;
1647:     PetscScalar *V;

1649:     PetscCall(PetscMalloc2(npoints, &diag, npoints, &subdiag));
1650:     PetscCall(PetscMalloc1(npoints * npoints, &V));
1651:     PetscCall(PetscDTJacobiMatrix_Internal(npoints, a, b, diag, subdiag));
1652:     for (i = 0; i < npoints - 1; i++) subdiag[i] = PetscSqrtReal(subdiag[i]);
1653:     PetscCall(PetscDTSymmetricTridiagonalEigensolve(npoints, diag, subdiag, x, V));
1654:     for (i = 0; i < npoints; i++) w[i] = PetscSqr(PetscRealPart(V[i * npoints])) * mu0;
1655:     PetscCall(PetscFree(V));
1656:     PetscCall(PetscFree2(diag, subdiag));
1657:   }
1658: #else
1659:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
1660: #endif
1661:   { /* As of March 2, 2020, The Sun Performance Library breaks the LAPACK contract for xstegr and xsteqr: the
1662:        eigenvalues are not guaranteed to be in ascending order.  So we heave a passive aggressive sigh and check that
1663:        the eigenvalues are sorted */
1664:     PetscBool sorted;

1666:     PetscCall(PetscSortedReal(npoints, x, &sorted));
1667:     if (!sorted) {
1668:       PetscInt  *order, i;
1669:       PetscReal *tmp;

1671:       PetscCall(PetscMalloc2(npoints, &order, npoints, &tmp));
1672:       for (i = 0; i < npoints; i++) order[i] = i;
1673:       PetscCall(PetscSortRealWithPermutation(npoints, x, order));
1674:       PetscCall(PetscArraycpy(tmp, x, npoints));
1675:       for (i = 0; i < npoints; i++) x[i] = tmp[order[i]];
1676:       PetscCall(PetscArraycpy(tmp, w, npoints));
1677:       for (i = 0; i < npoints; i++) w[i] = tmp[order[i]];
1678:       PetscCall(PetscFree2(order, tmp));
1679:     }
1680:   }
1681:   PetscFunctionReturn(PETSC_SUCCESS);
1682: }

1684: static PetscErrorCode PetscDTGaussJacobiQuadrature_Internal(PetscInt npoints, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
1685: {
1686:   PetscFunctionBegin;
1687:   PetscCheck(npoints >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be positive");
1688:   /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
1689:   PetscCheck(alpha > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "alpha must be > -1.");
1690:   PetscCheck(beta > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "beta must be > -1.");

1692:   if (newton) PetscCall(PetscDTGaussJacobiQuadrature_Newton_Internal(npoints, alpha, beta, x, w));
1693:   else PetscCall(PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(npoints, alpha, beta, x, w));
1694:   if (alpha == beta) { /* symmetrize */
1695:     PetscInt i;
1696:     for (i = 0; i < (npoints + 1) / 2; i++) {
1697:       PetscInt  j  = npoints - 1 - i;
1698:       PetscReal xi = x[i];
1699:       PetscReal xj = x[j];
1700:       PetscReal wi = w[i];
1701:       PetscReal wj = w[j];

1703:       x[i] = (xi - xj) / 2.;
1704:       x[j] = (xj - xi) / 2.;
1705:       w[i] = w[j] = (wi + wj) / 2.;
1706:     }
1707:   }
1708:   PetscFunctionReturn(PETSC_SUCCESS);
1709: }

1711: /*@
1712:   PetscDTGaussJacobiQuadrature - quadrature for the interval $[a, b]$ with the weight function
1713:   $(x-a)^\alpha (x-b)^\beta$.

1715:   Not Collective

1717:   Input Parameters:
1718: + npoints - the number of points in the quadrature rule
1719: . a       - the left endpoint of the interval
1720: . b       - the right endpoint of the interval
1721: . alpha   - the left exponent
1722: - beta    - the right exponent

1724:   Output Parameters:
1725: + x - array of length `npoints`, the locations of the quadrature points
1726: - w - array of length `npoints`, the weights of the quadrature points

1728:   Level: intermediate

1730:   Note:
1731:   This quadrature rule is exact for polynomials up to degree 2*`npoints` - 1.

1733: .seealso: `PetscDTGaussQuadrature()`
1734: @*/
1735: PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt npoints, PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
1736: {
1737:   PetscInt i;

1739:   PetscFunctionBegin;
1740:   PetscCall(PetscDTGaussJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal));
1741:   if (a != -1. || b != 1.) { /* shift */
1742:     for (i = 0; i < npoints; i++) {
1743:       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1744:       w[i] *= (b - a) / 2.;
1745:     }
1746:   }
1747:   PetscFunctionReturn(PETSC_SUCCESS);
1748: }

1750: static PetscErrorCode PetscDTGaussLobattoJacobiQuadrature_Internal(PetscInt npoints, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
1751: {
1752:   PetscFunctionBegin;
1753:   PetscCheck(npoints >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be positive");
1754:   /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
1755:   PetscCheck(alpha > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "alpha must be > -1.");
1756:   PetscCheck(beta > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "beta must be > -1.");

1758:   x[0]           = -1.;
1759:   x[npoints - 1] = 1.;
1760:   if (npoints > 2) PetscCall(PetscDTGaussJacobiQuadrature_Internal(npoints - 2, alpha + 1., beta + 1., &x[1], &w[1], newton));
1761:   for (PetscInt i = 1; i < npoints - 1; i++) w[i] /= (1. - x[i] * x[i]);
1762:   PetscCall(PetscDTGaussLobattoJacobiEndweights_Internal(npoints, alpha, beta, &w[0], &w[npoints - 1]));
1763:   PetscFunctionReturn(PETSC_SUCCESS);
1764: }

1766: /*@
1767:   PetscDTGaussLobattoJacobiQuadrature - quadrature for the interval $[a, b]$ with the weight function
1768:   $(x-a)^\alpha (x-b)^\beta$, with endpoints `a` and `b` included as quadrature points.

1770:   Not Collective

1772:   Input Parameters:
1773: + npoints - the number of points in the quadrature rule
1774: . a       - the left endpoint of the interval
1775: . b       - the right endpoint of the interval
1776: . alpha   - the left exponent
1777: - beta    - the right exponent

1779:   Output Parameters:
1780: + x - array of length `npoints`, the locations of the quadrature points
1781: - w - array of length `npoints`, the weights of the quadrature points

1783:   Level: intermediate

1785:   Note:
1786:   This quadrature rule is exact for polynomials up to degree 2*`npoints` - 3.

1788: .seealso: `PetscDTGaussJacobiQuadrature()`
1789: @*/
1790: PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt npoints, PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
1791: {
1792:   PetscFunctionBegin;
1793:   PetscCall(PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal));
1794:   if (a != -1. || b != 1.) { /* shift */
1795:     for (PetscInt i = 0; i < npoints; i++) {
1796:       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1797:       w[i] *= (b - a) / 2.;
1798:     }
1799:   }
1800:   PetscFunctionReturn(PETSC_SUCCESS);
1801: }

1803: /*@
1804:   PetscDTGaussQuadrature - create Gauss-Legendre quadrature

1806:   Not Collective

1808:   Input Parameters:
1809: + npoints - number of points
1810: . a       - left end of interval (often-1)
1811: - b       - right end of interval (often +1)

1813:   Output Parameters:
1814: + x - quadrature points
1815: - w - quadrature weights

1817:   Level: intermediate

1819:   Note:
1820:   See {cite}`golub1969calculation`

1822: .seealso: `PetscDTLegendreEval()`, `PetscDTGaussJacobiQuadrature()`
1823: @*/
1824: PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[])
1825: {
1826:   PetscFunctionBegin;
1827:   PetscCall(PetscDTGaussJacobiQuadrature_Internal(npoints, 0., 0., x, w, PetscDTGaussQuadratureNewton_Internal));
1828:   if (a != -1. || b != 1.) { /* shift */
1829:     for (PetscInt i = 0; i < npoints; i++) {
1830:       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1831:       w[i] *= (b - a) / 2.;
1832:     }
1833:   }
1834:   PetscFunctionReturn(PETSC_SUCCESS);
1835: }

1837: /*@
1838:   PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre
1839:   nodes of a given size on the domain $[-1,1]$

1841:   Not Collective

1843:   Input Parameters:
1844: + npoints - number of grid nodes
1845: - type    - `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` or `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON`

1847:   Output Parameters:
1848: + x - quadrature points, pass in an array of length `npoints`
1849: - w - quadrature weights, pass in an array of length `npoints`

1851:   Level: intermediate

1853:   Notes:
1854:   For n > 30  the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
1855:   close enough to the desired solution

1857:   These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes

1859:   See <https://epubs.siam.org/doi/abs/10.1137/110855442>  <https://epubs.siam.org/doi/abs/10.1137/120889873> for better ways to compute GLL nodes

1861: .seealso: `PetscDTGaussQuadrature()`, `PetscGaussLobattoLegendreCreateType`
1862: @*/
1863: PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints, PetscGaussLobattoLegendreCreateType type, PetscReal x[], PetscReal w[])
1864: {
1865:   PetscBool newton;

1867:   PetscFunctionBegin;
1868:   PetscCheck(npoints >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must provide at least 2 grid points per element");
1869:   newton = (PetscBool)(type == PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON);
1870:   PetscCall(PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, 0., 0., x, w, newton));
1871:   PetscFunctionReturn(PETSC_SUCCESS);
1872: }

1874: /*@
1875:   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature

1877:   Not Collective

1879:   Input Parameters:
1880: + dim     - The spatial dimension
1881: . Nc      - The number of components
1882: . npoints - number of points in one dimension
1883: . a       - left end of interval (often-1)
1884: - b       - right end of interval (often +1)

1886:   Output Parameter:
1887: . q - A `PetscQuadrature` object

1889:   Level: intermediate

1891: .seealso: `PetscDTGaussQuadrature()`, `PetscDTLegendreEval()`
1892: @*/
1893: PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1894: {
1895:   DMPolytopeType ct;
1896:   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints * PetscSqr(npoints) : PetscSqr(npoints) : npoints;
1897:   PetscReal     *x, *w, *xw, *ww;

1899:   PetscFunctionBegin;
1900:   PetscCall(PetscMalloc1(totpoints * dim, &x));
1901:   PetscCall(PetscMalloc1(totpoints * Nc, &w));
1902:   /* Set up the Golub-Welsch system */
1903:   switch (dim) {
1904:   case 0:
1905:     ct = DM_POLYTOPE_POINT;
1906:     PetscCall(PetscFree(x));
1907:     PetscCall(PetscFree(w));
1908:     PetscCall(PetscMalloc1(1, &x));
1909:     PetscCall(PetscMalloc1(Nc, &w));
1910:     totpoints = 1;
1911:     x[0]      = 0.0;
1912:     for (PetscInt c = 0; c < Nc; ++c) w[c] = 1.0;
1913:     break;
1914:   case 1:
1915:     ct = DM_POLYTOPE_SEGMENT;
1916:     PetscCall(PetscMalloc1(npoints, &ww));
1917:     PetscCall(PetscDTGaussQuadrature(npoints, a, b, x, ww));
1918:     for (PetscInt i = 0; i < npoints; ++i)
1919:       for (PetscInt c = 0; c < Nc; ++c) w[i * Nc + c] = ww[i];
1920:     PetscCall(PetscFree(ww));
1921:     break;
1922:   case 2:
1923:     ct = DM_POLYTOPE_QUADRILATERAL;
1924:     PetscCall(PetscMalloc2(npoints, &xw, npoints, &ww));
1925:     PetscCall(PetscDTGaussQuadrature(npoints, a, b, xw, ww));
1926:     for (PetscInt i = 0; i < npoints; ++i) {
1927:       for (PetscInt j = 0; j < npoints; ++j) {
1928:         x[(i * npoints + j) * dim + 0] = xw[i];
1929:         x[(i * npoints + j) * dim + 1] = xw[j];
1930:         for (PetscInt c = 0; c < Nc; ++c) w[(i * npoints + j) * Nc + c] = ww[i] * ww[j];
1931:       }
1932:     }
1933:     PetscCall(PetscFree2(xw, ww));
1934:     break;
1935:   case 3:
1936:     ct = DM_POLYTOPE_HEXAHEDRON;
1937:     PetscCall(PetscMalloc2(npoints, &xw, npoints, &ww));
1938:     PetscCall(PetscDTGaussQuadrature(npoints, a, b, xw, ww));
1939:     for (PetscInt i = 0; i < npoints; ++i) {
1940:       for (PetscInt j = 0; j < npoints; ++j) {
1941:         for (PetscInt k = 0; k < npoints; ++k) {
1942:           x[((i * npoints + j) * npoints + k) * dim + 0] = xw[i];
1943:           x[((i * npoints + j) * npoints + k) * dim + 1] = xw[j];
1944:           x[((i * npoints + j) * npoints + k) * dim + 2] = xw[k];
1945:           for (PetscInt c = 0; c < Nc; ++c) w[((i * npoints + j) * npoints + k) * Nc + c] = ww[i] * ww[j] * ww[k];
1946:         }
1947:       }
1948:     }
1949:     PetscCall(PetscFree2(xw, ww));
1950:     break;
1951:   default:
1952:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %" PetscInt_FMT, dim);
1953:   }
1954:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
1955:   PetscCall(PetscQuadratureSetCellType(*q, ct));
1956:   PetscCall(PetscQuadratureSetOrder(*q, 2 * npoints - 1));
1957:   PetscCall(PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w));
1958:   PetscCall(PetscObjectChangeTypeName((PetscObject)*q, "GaussTensor"));
1959:   PetscFunctionReturn(PETSC_SUCCESS);
1960: }

1962: /*@
1963:   PetscDTStroudConicalQuadrature - create Stroud conical quadrature for a simplex {cite}`karniadakis2005spectral`

1965:   Not Collective

1967:   Input Parameters:
1968: + dim     - The simplex dimension
1969: . Nc      - The number of components
1970: . npoints - The number of points in one dimension
1971: . a       - left end of interval (often-1)
1972: - b       - right end of interval (often +1)

1974:   Output Parameter:
1975: . q - A `PetscQuadrature` object

1977:   Level: intermediate

1979:   Note:
1980:   For `dim` == 1, this is Gauss-Legendre quadrature

1982: .seealso: `PetscDTGaussTensorQuadrature()`, `PetscDTGaussQuadrature()`
1983: @*/
1984: PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1985: {
1986:   DMPolytopeType ct;
1987:   PetscInt       totpoints;
1988:   PetscReal     *p1, *w1;
1989:   PetscReal     *x, *w;

1991:   PetscFunctionBegin;
1992:   PetscCheck(!(a != -1.0) && !(b != 1.0), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
1993:   switch (dim) {
1994:   case 0:
1995:     ct = DM_POLYTOPE_POINT;
1996:     break;
1997:   case 1:
1998:     ct = DM_POLYTOPE_SEGMENT;
1999:     break;
2000:   case 2:
2001:     ct = DM_POLYTOPE_TRIANGLE;
2002:     break;
2003:   case 3:
2004:     ct = DM_POLYTOPE_TETRAHEDRON;
2005:     break;
2006:   default:
2007:     ct = DM_POLYTOPE_UNKNOWN;
2008:   }
2009:   totpoints = 1;
2010:   for (PetscInt i = 0; i < dim; ++i) totpoints *= npoints;
2011:   PetscCall(PetscMalloc1(totpoints * dim, &x));
2012:   PetscCall(PetscMalloc1(totpoints * Nc, &w));
2013:   PetscCall(PetscMalloc2(npoints, &p1, npoints, &w1));
2014:   for (PetscInt i = 0; i < totpoints * Nc; ++i) w[i] = 1.;
2015:   for (PetscInt i = 0, totprev = 1, totrem = totpoints / npoints; i < dim; ++i) {
2016:     PetscReal mul;

2018:     mul = PetscPowReal(2., -i);
2019:     PetscCall(PetscDTGaussJacobiQuadrature(npoints, -1., 1., i, 0.0, p1, w1));
2020:     for (PetscInt pt = 0, l = 0; l < totprev; l++) {
2021:       for (PetscInt j = 0; j < npoints; j++) {
2022:         for (PetscInt m = 0; m < totrem; m++, pt++) {
2023:           for (PetscInt k = 0; k < i; k++) x[pt * dim + k] = (x[pt * dim + k] + 1.) * (1. - p1[j]) * 0.5 - 1.;
2024:           x[pt * dim + i] = p1[j];
2025:           for (PetscInt c = 0; c < Nc; c++) w[pt * Nc + c] *= mul * w1[j];
2026:         }
2027:       }
2028:     }
2029:     totprev *= npoints;
2030:     totrem /= npoints;
2031:   }
2032:   PetscCall(PetscFree2(p1, w1));
2033:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
2034:   PetscCall(PetscQuadratureSetCellType(*q, ct));
2035:   PetscCall(PetscQuadratureSetOrder(*q, 2 * npoints - 1));
2036:   PetscCall(PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w));
2037:   PetscCall(PetscObjectChangeTypeName((PetscObject)*q, "StroudConical"));
2038:   PetscFunctionReturn(PETSC_SUCCESS);
2039: }

2041: static PetscBool MinSymTriQuadCite       = PETSC_FALSE;
2042: const char       MinSymTriQuadCitation[] = "@article{WitherdenVincent2015,\n"
2043:                                            "  title = {On the identification of symmetric quadrature rules for finite element methods},\n"
2044:                                            "  journal = {Computers & Mathematics with Applications},\n"
2045:                                            "  volume = {69},\n"
2046:                                            "  number = {10},\n"
2047:                                            "  pages = {1232-1241},\n"
2048:                                            "  year = {2015},\n"
2049:                                            "  issn = {0898-1221},\n"
2050:                                            "  doi = {10.1016/j.camwa.2015.03.017},\n"
2051:                                            "  url = {https://www.sciencedirect.com/science/article/pii/S0898122115001224},\n"
2052:                                            "  author = {F.D. Witherden and P.E. Vincent},\n"
2053:                                            "}\n";

2055: #include "petscdttriquadrules.h"

2057: static PetscBool MinSymTetQuadCite       = PETSC_FALSE;
2058: const char       MinSymTetQuadCitation[] = "@article{JaskowiecSukumar2021\n"
2059:                                            "  author = {Jaskowiec, Jan and Sukumar, N.},\n"
2060:                                            "  title = {High-order symmetric cubature rules for tetrahedra and pyramids},\n"
2061:                                            "  journal = {International Journal for Numerical Methods in Engineering},\n"
2062:                                            "  volume = {122},\n"
2063:                                            "  number = {1},\n"
2064:                                            "  pages = {148-171},\n"
2065:                                            "  doi = {10.1002/nme.6528},\n"
2066:                                            "  url = {https://onlinelibrary.wiley.com/doi/abs/10.1002/nme.6528},\n"
2067:                                            "  eprint = {https://onlinelibrary.wiley.com/doi/pdf/10.1002/nme.6528},\n"
2068:                                            "  year = {2021}\n"
2069:                                            "}\n";

2071: #include "petscdttetquadrules.h"

2073: static PetscBool DiagSymTriQuadCite       = PETSC_FALSE;
2074: const char       DiagSymTriQuadCitation[] = "@article{KongMulderVeldhuizen1999,\n"
2075:                                             "  title = {Higher-order triangular and tetrahedral finite elements with mass lumping for solving the wave equation},\n"
2076:                                             "  journal = {Journal of Engineering Mathematics},\n"
2077:                                             "  volume = {35},\n"
2078:                                             "  number = {4},\n"
2079:                                             "  pages = {405--426},\n"
2080:                                             "  year = {1999},\n"
2081:                                             "  doi = {10.1023/A:1004420829610},\n"
2082:                                             "  url = {https://link.springer.com/article/10.1023/A:1004420829610},\n"
2083:                                             "  author = {MJS Chin-Joe-Kong and Wim A Mulder and Marinus Van Veldhuizen},\n"
2084:                                             "}\n";

2086: #include "petscdttridiagquadrules.h"

2088: // https://en.wikipedia.org/wiki/Partition_(number_theory)
2089: static PetscErrorCode PetscDTPartitionNumber(PetscInt n, PetscInt *p)
2090: {
2091:   // sequence A000041 in the OEIS
2092:   const PetscInt partition[]   = {1, 1, 2, 3, 5, 7, 11, 15, 22, 30, 42, 56, 77, 101, 135, 176, 231, 297, 385, 490, 627, 792, 1002, 1255, 1575, 1958, 2436, 3010, 3718, 4565, 5604};
2093:   PetscInt       tabulated_max = PETSC_STATIC_ARRAY_LENGTH(partition) - 1;

2095:   PetscFunctionBegin;
2096:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Partition number not defined for negative number %" PetscInt_FMT, n);
2097:   // not implementing the pentagonal number recurrence, we don't need partition numbers for n that high
2098:   PetscCheck(n <= tabulated_max, PETSC_COMM_SELF, PETSC_ERR_SUP, "Partition numbers only tabulated up to %" PetscInt_FMT ", not computed for %" PetscInt_FMT, tabulated_max, n);
2099:   *p = partition[n];
2100:   PetscFunctionReturn(PETSC_SUCCESS);
2101: }

2103: /*@
2104:   PetscDTSimplexQuadrature - Create a quadrature rule for a simplex that exactly integrates polynomials up to a given degree.

2106:   Not Collective

2108:   Input Parameters:
2109: + dim    - The spatial dimension of the simplex (1 = segment, 2 = triangle, 3 = tetrahedron)
2110: . degree - The largest polynomial degree that is required to be integrated exactly
2111: - type   - `PetscDTSimplexQuadratureType` indicating the type of quadrature rule

2113:   Output Parameter:
2114: . quad - A `PetscQuadrature` object for integration over the biunit simplex
2115:             (defined by the bounds $x_i >= -1$ and $\sum_i x_i <= 2 - d$) that is exact for
2116:             polynomials up to the given degree

2118:   Level: intermediate

2120: .seealso: `PetscDTSimplexQuadratureType`, `PetscDTGaussQuadrature()`, `PetscDTStroudCononicalQuadrature()`, `PetscQuadrature`
2121: @*/
2122: PetscErrorCode PetscDTSimplexQuadrature(PetscInt dim, PetscInt degree, PetscDTSimplexQuadratureType type, PetscQuadrature *quad)
2123: {
2124:   PetscDTSimplexQuadratureType orig_type = type;

2126:   PetscFunctionBegin;
2127:   PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative dimension %" PetscInt_FMT, dim);
2128:   PetscCheck(degree >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %" PetscInt_FMT, degree);
2129:   if (type == PETSCDTSIMPLEXQUAD_DEFAULT) type = PETSCDTSIMPLEXQUAD_MINSYM;
2130:   if (type == PETSCDTSIMPLEXQUAD_CONIC || dim < 2) {
2131:     PetscInt points_per_dim = (degree + 2) / 2; // ceil((degree + 1) / 2);
2132:     PetscCall(PetscDTStroudConicalQuadrature(dim, 1, points_per_dim, -1, 1, quad));
2133:   } else {
2134:     DMPolytopeType    ct;
2135:     PetscInt          n    = dim + 1;
2136:     PetscInt          fact = 1;
2137:     PetscInt         *part, *perm;
2138:     PetscInt          p = 0;
2139:     PetscInt          max_degree;
2140:     const PetscInt   *nodes_per_type     = NULL;
2141:     const PetscInt   *all_num_full_nodes = NULL;
2142:     const PetscReal **weights_list       = NULL;
2143:     const PetscReal **compact_nodes_list = NULL;
2144:     const char       *citation           = NULL;
2145:     PetscBool        *cited              = NULL;

2147:     switch (dim) {
2148:     case 0:
2149:       ct = DM_POLYTOPE_POINT;
2150:       break;
2151:     case 1:
2152:       ct = DM_POLYTOPE_SEGMENT;
2153:       break;
2154:     case 2:
2155:       ct = DM_POLYTOPE_TRIANGLE;
2156:       break;
2157:     case 3:
2158:       ct = DM_POLYTOPE_TETRAHEDRON;
2159:       break;
2160:     default:
2161:       ct = DM_POLYTOPE_UNKNOWN;
2162:     }
2163:     if (type == PETSCDTSIMPLEXQUAD_MINSYM) {
2164:       switch (dim) {
2165:       case 2:
2166:         cited              = &MinSymTriQuadCite;
2167:         citation           = MinSymTriQuadCitation;
2168:         max_degree         = PetscDTWVTriQuad_max_degree;
2169:         nodes_per_type     = PetscDTWVTriQuad_num_orbits;
2170:         all_num_full_nodes = PetscDTWVTriQuad_num_nodes;
2171:         weights_list       = PetscDTWVTriQuad_weights;
2172:         compact_nodes_list = PetscDTWVTriQuad_orbits;
2173:         break;
2174:       case 3:
2175:         cited              = &MinSymTetQuadCite;
2176:         citation           = MinSymTetQuadCitation;
2177:         max_degree         = PetscDTJSTetQuad_max_degree;
2178:         nodes_per_type     = PetscDTJSTetQuad_num_orbits;
2179:         all_num_full_nodes = PetscDTJSTetQuad_num_nodes;
2180:         weights_list       = PetscDTJSTetQuad_weights;
2181:         compact_nodes_list = PetscDTJSTetQuad_orbits;
2182:         break;
2183:       default:
2184:         max_degree = -1;
2185:         break;
2186:       }
2187:     } else {
2188:       switch (dim) {
2189:       case 2:
2190:         cited              = &DiagSymTriQuadCite;
2191:         citation           = DiagSymTriQuadCitation;
2192:         max_degree         = PetscDTKMVTriQuad_max_degree;
2193:         nodes_per_type     = PetscDTKMVTriQuad_num_orbits;
2194:         all_num_full_nodes = PetscDTKMVTriQuad_num_nodes;
2195:         weights_list       = PetscDTKMVTriQuad_weights;
2196:         compact_nodes_list = PetscDTKMVTriQuad_orbits;
2197:         break;
2198:       default:
2199:         max_degree = -1;
2200:         break;
2201:       }
2202:     }

2204:     if (degree > max_degree) {
2205:       PetscCheck(orig_type == PETSCDTSIMPLEXQUAD_DEFAULT, PETSC_COMM_SELF, PETSC_ERR_SUP, "%s symmetric quadrature for dim %" PetscInt_FMT ", degree %" PetscInt_FMT " unsupported", orig_type == PETSCDTSIMPLEXQUAD_MINSYM ? "Minimal" : "Diagonal", dim, degree);
2206:       // fall back to conic
2207:       PetscCall(PetscDTSimplexQuadrature(dim, degree, PETSCDTSIMPLEXQUAD_CONIC, quad));
2208:       PetscFunctionReturn(PETSC_SUCCESS);
2209:     }

2211:     PetscCall(PetscCitationsRegister(citation, cited));

2213:     PetscCall(PetscDTPartitionNumber(n, &p));
2214:     for (PetscInt d = 2; d <= n; d++) fact *= d;

2216:     PetscInt         num_full_nodes      = all_num_full_nodes[degree];
2217:     const PetscReal *all_compact_nodes   = compact_nodes_list[degree];
2218:     const PetscReal *all_compact_weights = weights_list[degree];
2219:     nodes_per_type                       = &nodes_per_type[p * degree];

2221:     PetscReal      *points;
2222:     PetscReal      *counts;
2223:     PetscReal      *weights;
2224:     PetscReal      *bary_to_biunit; // row-major transformation of barycentric coordinate to biunit
2225:     PetscQuadrature q;

2227:     // compute the transformation
2228:     PetscCall(PetscMalloc1(n * dim, &bary_to_biunit));
2229:     for (PetscInt d = 0; d < dim; d++) {
2230:       for (PetscInt b = 0; b < n; b++) bary_to_biunit[d * n + b] = (d == b) ? 1.0 : -1.0;
2231:     }

2233:     PetscCall(PetscMalloc3(n, &part, n, &perm, n, &counts));
2234:     PetscCall(PetscCalloc1(num_full_nodes * dim, &points));
2235:     PetscCall(PetscMalloc1(num_full_nodes, &weights));

2237:     // (0, 0, ...) is the first partition lexicographically
2238:     PetscCall(PetscArrayzero(part, n));
2239:     PetscCall(PetscArrayzero(counts, n));
2240:     counts[0] = n;

2242:     // for each partition
2243:     for (PetscInt s = 0, node_offset = 0; s < p; s++) {
2244:       PetscInt num_compact_coords = part[n - 1] + 1;

2246:       const PetscReal *compact_nodes   = all_compact_nodes;
2247:       const PetscReal *compact_weights = all_compact_weights;
2248:       all_compact_nodes += num_compact_coords * nodes_per_type[s];
2249:       all_compact_weights += nodes_per_type[s];

2251:       // for every permutation of the vertices
2252:       for (PetscInt f = 0; f < fact; f++) {
2253:         PetscCall(PetscDTEnumPerm(n, f, perm, NULL));

2255:         // check if it is a valid permutation
2256:         PetscInt digit;
2257:         for (digit = 1; digit < n; digit++) {
2258:           // skip permutations that would duplicate a node because it has a smaller symmetry group
2259:           if (part[digit - 1] == part[digit] && perm[digit - 1] > perm[digit]) break;
2260:         }
2261:         if (digit < n) continue;

2263:         // create full nodes from this permutation of the compact nodes
2264:         PetscReal *full_nodes   = &points[node_offset * dim];
2265:         PetscReal *full_weights = &weights[node_offset];

2267:         PetscCall(PetscArraycpy(full_weights, compact_weights, nodes_per_type[s]));
2268:         for (PetscInt b = 0; b < n; b++) {
2269:           for (PetscInt d = 0; d < dim; d++) {
2270:             for (PetscInt node = 0; node < nodes_per_type[s]; node++) full_nodes[node * dim + d] += bary_to_biunit[d * n + perm[b]] * compact_nodes[node * num_compact_coords + part[b]];
2271:           }
2272:         }
2273:         node_offset += nodes_per_type[s];
2274:       }

2276:       if (s < p - 1) { // Generate the next partition
2277:         /* A partition is described by the number of coordinates that are in
2278:          * each set of duplicates (counts) and redundantly by mapping each
2279:          * index to its set of duplicates (part)
2280:          *
2281:          * Counts should always be in nonincreasing order
2282:          *
2283:          * We want to generate the partitions lexically by part, which means
2284:          * finding the last index where count > 1 and reducing by 1.
2285:          *
2286:          * For the new counts beyond that index, we eagerly assign the remaining
2287:          * capacity of the partition to smaller indices (ensures lexical ordering),
2288:          * while respecting the nonincreasing invariant of the counts
2289:          */
2290:         PetscInt last_digit            = part[n - 1];
2291:         PetscInt last_digit_with_extra = last_digit;
2292:         while (counts[last_digit_with_extra] == 1) last_digit_with_extra--;
2293:         PetscInt limit               = --counts[last_digit_with_extra];
2294:         PetscInt total_to_distribute = last_digit - last_digit_with_extra + 1;
2295:         for (PetscInt digit = last_digit_with_extra + 1; digit < n; digit++) {
2296:           counts[digit] = PetscMin(limit, total_to_distribute);
2297:           total_to_distribute -= PetscMin(limit, total_to_distribute);
2298:         }
2299:         for (PetscInt digit = 0, offset = 0; digit < n; digit++) {
2300:           PetscInt count = counts[digit];
2301:           for (PetscInt c = 0; c < count; c++) part[offset++] = digit;
2302:         }
2303:       }
2304:       PetscCheck(node_offset <= num_full_nodes, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node offset %" PetscInt_FMT " > %" PetscInt_FMT " number of nodes", node_offset, num_full_nodes);
2305:     }
2306:     PetscCall(PetscFree3(part, perm, counts));
2307:     PetscCall(PetscFree(bary_to_biunit));
2308:     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &q));
2309:     PetscCall(PetscQuadratureSetCellType(q, ct));
2310:     PetscCall(PetscQuadratureSetOrder(q, degree));
2311:     PetscCall(PetscQuadratureSetData(q, dim, 1, num_full_nodes, points, weights));
2312:     *quad = q;
2313:   }
2314:   PetscFunctionReturn(PETSC_SUCCESS);
2315: }

2317: /*@
2318:   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell

2320:   Not Collective

2322:   Input Parameters:
2323: + dim   - The cell dimension
2324: . level - The number of points in one dimension, $2^l$
2325: . a     - left end of interval (often-1)
2326: - b     - right end of interval (often +1)

2328:   Output Parameter:
2329: . q - A `PetscQuadrature` object

2331:   Level: intermediate

2333: .seealso: `PetscDTGaussTensorQuadrature()`, `PetscQuadrature`
2334: @*/
2335: PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
2336: {
2337:   DMPolytopeType  ct;
2338:   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
2339:   const PetscReal alpha = (b - a) / 2.;              /* Half-width of the integration interval */
2340:   const PetscReal beta  = (b + a) / 2.;              /* Center of the integration interval */
2341:   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
2342:   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
2343:   PetscReal       wk = 0.5 * PETSC_PI;               /* Quadrature weight at x_k */
2344:   PetscReal      *x, *w;
2345:   PetscInt        K, k, npoints;

2347:   PetscFunctionBegin;
2348:   PetscCheck(dim <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not yet implemented", dim);
2349:   PetscCheck(level, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
2350:   switch (dim) {
2351:   case 0:
2352:     ct = DM_POLYTOPE_POINT;
2353:     break;
2354:   case 1:
2355:     ct = DM_POLYTOPE_SEGMENT;
2356:     break;
2357:   case 2:
2358:     ct = DM_POLYTOPE_QUADRILATERAL;
2359:     break;
2360:   case 3:
2361:     ct = DM_POLYTOPE_HEXAHEDRON;
2362:     break;
2363:   default:
2364:     ct = DM_POLYTOPE_UNKNOWN;
2365:   }
2366:   /* Find K such that the weights are < 32 digits of precision */
2367:   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2 * p; ++K) wk = 0.5 * h * PETSC_PI * PetscCoshReal(K * h) / PetscSqr(PetscCoshReal(0.5 * PETSC_PI * PetscSinhReal(K * h)));
2368:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
2369:   PetscCall(PetscQuadratureSetCellType(*q, ct));
2370:   PetscCall(PetscQuadratureSetOrder(*q, 2 * K + 1));
2371:   npoints = 2 * K - 1;
2372:   PetscCall(PetscMalloc1(npoints * dim, &x));
2373:   PetscCall(PetscMalloc1(npoints, &w));
2374:   /* Center term */
2375:   x[0] = beta;
2376:   w[0] = 0.5 * alpha * PETSC_PI;
2377:   for (k = 1; k < K; ++k) {
2378:     wk           = 0.5 * alpha * h * PETSC_PI * PetscCoshReal(k * h) / PetscSqr(PetscCoshReal(0.5 * PETSC_PI * PetscSinhReal(k * h)));
2379:     xk           = PetscTanhReal(0.5 * PETSC_PI * PetscSinhReal(k * h));
2380:     x[2 * k - 1] = -alpha * xk + beta;
2381:     w[2 * k - 1] = wk;
2382:     x[2 * k + 0] = alpha * xk + beta;
2383:     w[2 * k + 0] = wk;
2384:   }
2385:   PetscCall(PetscQuadratureSetData(*q, dim, 1, npoints, x, w));
2386:   PetscFunctionReturn(PETSC_SUCCESS);
2387: }

2389: PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(const PetscReal[], void *, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscCtx ctx, PetscReal *sol)
2390: {
2391:   const PetscInt  p     = 16;           /* Digits of precision in the evaluation */
2392:   const PetscReal alpha = (b - a) / 2.; /* Half-width of the integration interval */
2393:   const PetscReal beta  = (b + a) / 2.; /* Center of the integration interval */
2394:   PetscReal       h     = 1.0;          /* Step size, length between x_k */
2395:   PetscInt        l     = 0;            /* Level of refinement, h = 2^{-l} */
2396:   PetscReal       osum  = 0.0;          /* Integral on last level */
2397:   PetscReal       psum  = 0.0;          /* Integral on the level before the last level */
2398:   PetscReal       sum;                  /* Integral on current level */
2399:   PetscReal       yk;                   /* Quadrature point 1 - x_k on reference domain [-1, 1] */
2400:   PetscReal       lx, rx;               /* Quadrature points to the left and right of 0 on the real domain [a, b] */
2401:   PetscReal       wk;                   /* Quadrature weight at x_k */
2402:   PetscReal       lval, rval;           /* Terms in the quadature sum to the left and right of 0 */
2403:   PetscInt        d;                    /* Digits of precision in the integral */

2405:   PetscFunctionBegin;
2406:   PetscCheck(digits > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
2407:   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
2408:   /* Center term */
2409:   func(&beta, ctx, &lval);
2410:   sum = 0.5 * alpha * PETSC_PI * lval;
2411:   /* */
2412:   do {
2413:     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
2414:     PetscInt  k = 1;

2416:     ++l;
2417:     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %" PetscInt_FMT " sum: %15.15f\n", l, sum); */
2418:     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
2419:     psum = osum;
2420:     osum = sum;
2421:     h *= 0.5;
2422:     sum *= 0.5;
2423:     do {
2424:       wk = 0.5 * h * PETSC_PI * PetscCoshReal(k * h) / PetscSqr(PetscCoshReal(0.5 * PETSC_PI * PetscSinhReal(k * h)));
2425:       yk = 1.0 / (PetscExpReal(0.5 * PETSC_PI * PetscSinhReal(k * h)) * PetscCoshReal(0.5 * PETSC_PI * PetscSinhReal(k * h)));
2426:       lx = -alpha * (1.0 - yk) + beta;
2427:       rx = alpha * (1.0 - yk) + beta;
2428:       func(&lx, ctx, &lval);
2429:       func(&rx, ctx, &rval);
2430:       lterm   = alpha * wk * lval;
2431:       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
2432:       sum += lterm;
2433:       rterm   = alpha * wk * rval;
2434:       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
2435:       sum += rterm;
2436:       ++k;
2437:       /* Only need to evaluate every other point on refined levels */
2438:       if (l != 1) ++k;
2439:     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */

2441:     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
2442:     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
2443:     d3 = PetscLog10Real(maxTerm) - p;
2444:     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
2445:     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
2446:     d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1) / d2, 2 * d1), d3), d4)));
2447:   } while (d < digits && l < 12);
2448:   *sol = sum;
2449:   PetscCall(PetscFPTrapPop());
2450:   PetscFunctionReturn(PETSC_SUCCESS);
2451: }

2453: #if defined(PETSC_HAVE_MPFR)
2454: PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(const PetscReal[], void *, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscCtx ctx, PetscReal *sol)
2455: {
2456:   const PetscInt safetyFactor = 2; /* Calculate abscissa until 2*p digits */
2457:   PetscInt       l            = 0; /* Level of refinement, h = 2^{-l} */
2458:   mpfr_t         alpha;            /* Half-width of the integration interval */
2459:   mpfr_t         beta;             /* Center of the integration interval */
2460:   mpfr_t         h;                /* Step size, length between x_k */
2461:   mpfr_t         osum;             /* Integral on last level */
2462:   mpfr_t         psum;             /* Integral on the level before the last level */
2463:   mpfr_t         sum;              /* Integral on current level */
2464:   mpfr_t         yk;               /* Quadrature point 1 - x_k on reference domain [-1, 1] */
2465:   mpfr_t         lx, rx;           /* Quadrature points to the left and right of 0 on the real domain [a, b] */
2466:   mpfr_t         wk;               /* Quadrature weight at x_k */
2467:   PetscReal      lval, rval, rtmp; /* Terms in the quadature sum to the left and right of 0 */
2468:   PetscInt       d;                /* Digits of precision in the integral */
2469:   mpfr_t         pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;

2471:   PetscFunctionBegin;
2472:   PetscCheck(digits > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
2473:   /* Create high precision storage */
2474:   mpfr_inits2(PetscCeilReal(safetyFactor * digits * PetscLogReal(10.) / PetscLogReal(2.)), alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
2475:   /* Initialization */
2476:   mpfr_set_d(alpha, 0.5 * (b - a), MPFR_RNDN);
2477:   mpfr_set_d(beta, 0.5 * (b + a), MPFR_RNDN);
2478:   mpfr_set_d(osum, 0.0, MPFR_RNDN);
2479:   mpfr_set_d(psum, 0.0, MPFR_RNDN);
2480:   mpfr_set_d(h, 1.0, MPFR_RNDN);
2481:   mpfr_const_pi(pi2, MPFR_RNDN);
2482:   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
2483:   /* Center term */
2484:   rtmp = 0.5 * (b + a);
2485:   func(&rtmp, ctx, &lval);
2486:   mpfr_set(sum, pi2, MPFR_RNDN);
2487:   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
2488:   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
2489:   /* */
2490:   do {
2491:     PetscReal d1, d2, d3, d4;
2492:     PetscInt  k = 1;

2494:     ++l;
2495:     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
2496:     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %" PetscInt_FMT " sum: %15.15f\n", l, sum); */
2497:     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
2498:     mpfr_set(psum, osum, MPFR_RNDN);
2499:     mpfr_set(osum, sum, MPFR_RNDN);
2500:     mpfr_mul_d(h, h, 0.5, MPFR_RNDN);
2501:     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
2502:     do {
2503:       mpfr_set_si(kh, k, MPFR_RNDN);
2504:       mpfr_mul(kh, kh, h, MPFR_RNDN);
2505:       /* Weight */
2506:       mpfr_set(wk, h, MPFR_RNDN);
2507:       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
2508:       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
2509:       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
2510:       mpfr_cosh(tmp, msinh, MPFR_RNDN);
2511:       mpfr_sqr(tmp, tmp, MPFR_RNDN);
2512:       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
2513:       mpfr_div(wk, wk, tmp, MPFR_RNDN);
2514:       /* Abscissa */
2515:       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
2516:       mpfr_cosh(tmp, msinh, MPFR_RNDN);
2517:       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
2518:       mpfr_exp(tmp, msinh, MPFR_RNDN);
2519:       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
2520:       /* Quadrature points */
2521:       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
2522:       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
2523:       mpfr_add(lx, lx, beta, MPFR_RNDU);
2524:       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
2525:       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
2526:       mpfr_add(rx, rx, beta, MPFR_RNDD);
2527:       /* Evaluation */
2528:       rtmp = mpfr_get_d(lx, MPFR_RNDU);
2529:       func(&rtmp, ctx, &lval);
2530:       rtmp = mpfr_get_d(rx, MPFR_RNDD);
2531:       func(&rtmp, ctx, &rval);
2532:       /* Update */
2533:       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
2534:       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
2535:       mpfr_add(sum, sum, tmp, MPFR_RNDN);
2536:       mpfr_abs(tmp, tmp, MPFR_RNDN);
2537:       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
2538:       mpfr_set(curTerm, tmp, MPFR_RNDN);
2539:       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
2540:       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
2541:       mpfr_add(sum, sum, tmp, MPFR_RNDN);
2542:       mpfr_abs(tmp, tmp, MPFR_RNDN);
2543:       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
2544:       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
2545:       ++k;
2546:       /* Only need to evaluate every other point on refined levels */
2547:       if (l != 1) ++k;
2548:       mpfr_log10(tmp, wk, MPFR_RNDN);
2549:       mpfr_abs(tmp, tmp, MPFR_RNDN);
2550:     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor * digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
2551:     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
2552:     mpfr_abs(tmp, tmp, MPFR_RNDN);
2553:     mpfr_log10(tmp, tmp, MPFR_RNDN);
2554:     d1 = mpfr_get_d(tmp, MPFR_RNDN);
2555:     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
2556:     mpfr_abs(tmp, tmp, MPFR_RNDN);
2557:     mpfr_log10(tmp, tmp, MPFR_RNDN);
2558:     d2 = mpfr_get_d(tmp, MPFR_RNDN);
2559:     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
2560:     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
2561:     mpfr_log10(tmp, curTerm, MPFR_RNDN);
2562:     d4 = mpfr_get_d(tmp, MPFR_RNDN);
2563:     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1) / d2, 2 * d1), d3), d4)));
2564:   } while (d < digits && l < 8);
2565:   *sol = mpfr_get_d(sum, MPFR_RNDN);
2566:   /* Cleanup */
2567:   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
2568:   PetscFunctionReturn(PETSC_SUCCESS);
2569: }
2570: #else

2572: PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(const PetscReal[], void *, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscCtx ctx, PetscReal *sol)
2573: {
2574:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
2575: }
2576: #endif

2578: /*@
2579:   PetscDTTensorQuadratureCreate - create the tensor product quadrature from two lower-dimensional quadratures

2581:   Not Collective

2583:   Input Parameters:
2584: + q1 - The first quadrature
2585: - q2 - The second quadrature

2587:   Output Parameter:
2588: . q - A `PetscQuadrature` object

2590:   Level: intermediate

2592: .seealso: `PetscQuadrature`, `PetscDTGaussTensorQuadrature()`
2593: @*/
2594: PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature q1, PetscQuadrature q2, PetscQuadrature *q)
2595: {
2596:   DMPolytopeType   ct1, ct2, ct;
2597:   const PetscReal *x1, *w1, *x2, *w2;
2598:   PetscReal       *x, *w;
2599:   PetscInt         dim1, Nc1, Np1, order1, qa, d1;
2600:   PetscInt         dim2, Nc2, Np2, order2, qb, d2;
2601:   PetscInt         dim, Nc, Np, order, qc, d;

2603:   PetscFunctionBegin;
2606:   PetscAssertPointer(q, 3);

2608:   PetscCall(PetscQuadratureGetOrder(q1, &order1));
2609:   PetscCall(PetscQuadratureGetOrder(q2, &order2));
2610:   PetscCheck(order1 == order2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Order1 %" PetscInt_FMT " != %" PetscInt_FMT " Order2", order1, order2);
2611:   PetscCall(PetscQuadratureGetData(q1, &dim1, &Nc1, &Np1, &x1, &w1));
2612:   PetscCall(PetscQuadratureGetCellType(q1, &ct1));
2613:   PetscCall(PetscQuadratureGetData(q2, &dim2, &Nc2, &Np2, &x2, &w2));
2614:   PetscCall(PetscQuadratureGetCellType(q2, &ct2));
2615:   PetscCheck(Nc1 == Nc2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "NumComp1 %" PetscInt_FMT " != %" PetscInt_FMT " NumComp2", Nc1, Nc2);

2617:   switch (ct1) {
2618:   case DM_POLYTOPE_POINT:
2619:     ct = ct2;
2620:     break;
2621:   case DM_POLYTOPE_SEGMENT:
2622:     switch (ct2) {
2623:     case DM_POLYTOPE_POINT:
2624:       ct = ct1;
2625:       break;
2626:     case DM_POLYTOPE_SEGMENT:
2627:       ct = DM_POLYTOPE_QUADRILATERAL;
2628:       break;
2629:     case DM_POLYTOPE_TRIANGLE:
2630:       ct = DM_POLYTOPE_TRI_PRISM;
2631:       break;
2632:     case DM_POLYTOPE_QUADRILATERAL:
2633:       ct = DM_POLYTOPE_HEXAHEDRON;
2634:       break;
2635:     case DM_POLYTOPE_TETRAHEDRON:
2636:       ct = DM_POLYTOPE_UNKNOWN;
2637:       break;
2638:     case DM_POLYTOPE_HEXAHEDRON:
2639:       ct = DM_POLYTOPE_UNKNOWN;
2640:       break;
2641:     default:
2642:       ct = DM_POLYTOPE_UNKNOWN;
2643:     }
2644:     break;
2645:   case DM_POLYTOPE_TRIANGLE:
2646:     switch (ct2) {
2647:     case DM_POLYTOPE_POINT:
2648:       ct = ct1;
2649:       break;
2650:     case DM_POLYTOPE_SEGMENT:
2651:       ct = DM_POLYTOPE_TRI_PRISM;
2652:       break;
2653:     case DM_POLYTOPE_TRIANGLE:
2654:       ct = DM_POLYTOPE_UNKNOWN;
2655:       break;
2656:     case DM_POLYTOPE_QUADRILATERAL:
2657:       ct = DM_POLYTOPE_UNKNOWN;
2658:       break;
2659:     case DM_POLYTOPE_TETRAHEDRON:
2660:       ct = DM_POLYTOPE_UNKNOWN;
2661:       break;
2662:     case DM_POLYTOPE_HEXAHEDRON:
2663:       ct = DM_POLYTOPE_UNKNOWN;
2664:       break;
2665:     default:
2666:       ct = DM_POLYTOPE_UNKNOWN;
2667:     }
2668:     break;
2669:   case DM_POLYTOPE_QUADRILATERAL:
2670:     switch (ct2) {
2671:     case DM_POLYTOPE_POINT:
2672:       ct = ct1;
2673:       break;
2674:     case DM_POLYTOPE_SEGMENT:
2675:       ct = DM_POLYTOPE_HEXAHEDRON;
2676:       break;
2677:     case DM_POLYTOPE_TRIANGLE:
2678:       ct = DM_POLYTOPE_UNKNOWN;
2679:       break;
2680:     case DM_POLYTOPE_QUADRILATERAL:
2681:       ct = DM_POLYTOPE_UNKNOWN;
2682:       break;
2683:     case DM_POLYTOPE_TETRAHEDRON:
2684:       ct = DM_POLYTOPE_UNKNOWN;
2685:       break;
2686:     case DM_POLYTOPE_HEXAHEDRON:
2687:       ct = DM_POLYTOPE_UNKNOWN;
2688:       break;
2689:     default:
2690:       ct = DM_POLYTOPE_UNKNOWN;
2691:     }
2692:     break;
2693:   case DM_POLYTOPE_TETRAHEDRON:
2694:     switch (ct2) {
2695:     case DM_POLYTOPE_POINT:
2696:       ct = ct1;
2697:       break;
2698:     case DM_POLYTOPE_SEGMENT:
2699:       ct = DM_POLYTOPE_UNKNOWN;
2700:       break;
2701:     case DM_POLYTOPE_TRIANGLE:
2702:       ct = DM_POLYTOPE_UNKNOWN;
2703:       break;
2704:     case DM_POLYTOPE_QUADRILATERAL:
2705:       ct = DM_POLYTOPE_UNKNOWN;
2706:       break;
2707:     case DM_POLYTOPE_TETRAHEDRON:
2708:       ct = DM_POLYTOPE_UNKNOWN;
2709:       break;
2710:     case DM_POLYTOPE_HEXAHEDRON:
2711:       ct = DM_POLYTOPE_UNKNOWN;
2712:       break;
2713:     default:
2714:       ct = DM_POLYTOPE_UNKNOWN;
2715:     }
2716:     break;
2717:   case DM_POLYTOPE_HEXAHEDRON:
2718:     switch (ct2) {
2719:     case DM_POLYTOPE_POINT:
2720:       ct = ct1;
2721:       break;
2722:     case DM_POLYTOPE_SEGMENT:
2723:       ct = DM_POLYTOPE_UNKNOWN;
2724:       break;
2725:     case DM_POLYTOPE_TRIANGLE:
2726:       ct = DM_POLYTOPE_UNKNOWN;
2727:       break;
2728:     case DM_POLYTOPE_QUADRILATERAL:
2729:       ct = DM_POLYTOPE_UNKNOWN;
2730:       break;
2731:     case DM_POLYTOPE_TETRAHEDRON:
2732:       ct = DM_POLYTOPE_UNKNOWN;
2733:       break;
2734:     case DM_POLYTOPE_HEXAHEDRON:
2735:       ct = DM_POLYTOPE_UNKNOWN;
2736:       break;
2737:     default:
2738:       ct = DM_POLYTOPE_UNKNOWN;
2739:     }
2740:     break;
2741:   default:
2742:     ct = DM_POLYTOPE_UNKNOWN;
2743:   }
2744:   dim   = dim1 + dim2;
2745:   Nc    = Nc1;
2746:   Np    = Np1 * Np2;
2747:   order = order1;
2748:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
2749:   PetscCall(PetscQuadratureSetCellType(*q, ct));
2750:   PetscCall(PetscQuadratureSetOrder(*q, order));
2751:   PetscCall(PetscMalloc1(Np * dim, &x));
2752:   PetscCall(PetscMalloc1(Np, &w));
2753:   for (qa = 0, qc = 0; qa < Np1; ++qa) {
2754:     for (qb = 0; qb < Np2; ++qb, ++qc) {
2755:       for (d1 = 0, d = 0; d1 < dim1; ++d1, ++d) x[qc * dim + d] = x1[qa * dim1 + d1];
2756:       for (d2 = 0; d2 < dim2; ++d2, ++d) x[qc * dim + d] = x2[qb * dim2 + d2];
2757:       w[qc] = w1[qa] * w2[qb];
2758:     }
2759:   }
2760:   PetscCall(PetscQuadratureSetData(*q, dim, Nc, Np, x, w));
2761:   PetscFunctionReturn(PETSC_SUCCESS);
2762: }

2764: /* Overwrites A. Can only handle full-rank problems with m>=n
2765:    A in column-major format
2766:    Ainv in row-major format
2767:    tau has length m
2768:    worksize must be >= max(1,n)
2769:  */
2770: static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m, PetscInt mstride, PetscInt n, PetscReal *A_in, PetscReal *Ainv_out, PetscScalar *tau, PetscInt worksize, PetscScalar *work)
2771: {
2772:   PetscBLASInt M, N, K, lda, ldb, ldwork, info;
2773:   PetscScalar *A, *Ainv, *R, *Q, Alpha;

2775:   PetscFunctionBegin;
2776: #if defined(PETSC_USE_COMPLEX)
2777:   {
2778:     PetscInt i, j;
2779:     PetscCall(PetscMalloc2(m * n, &A, m * n, &Ainv));
2780:     for (j = 0; j < n; j++) {
2781:       for (i = 0; i < m; i++) A[i + m * j] = A_in[i + mstride * j];
2782:     }
2783:     mstride = m;
2784:   }
2785: #else
2786:   A    = A_in;
2787:   Ainv = Ainv_out;
2788: #endif

2790:   PetscCall(PetscBLASIntCast(m, &M));
2791:   PetscCall(PetscBLASIntCast(n, &N));
2792:   PetscCall(PetscBLASIntCast(mstride, &lda));
2793:   PetscCall(PetscBLASIntCast(worksize, &ldwork));
2794:   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
2795:   PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&M, &N, A, &lda, tau, work, &ldwork, &info));
2796:   PetscCall(PetscFPTrapPop());
2797:   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "xGEQRF error");
2798:   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */

2800:   /* Extract an explicit representation of Q */
2801:   Q = Ainv;
2802:   PetscCall(PetscArraycpy(Q, A, mstride * n));
2803:   K = N; /* full rank */
2804:   PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&M, &N, &K, Q, &lda, tau, work, &ldwork, &info));
2805:   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "xORGQR/xUNGQR error");

2807:   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
2808:   Alpha = 1.0;
2809:   ldb   = lda;
2810:   PetscCallBLAS("BLAStrsm", BLAStrsm_("Right", "Upper", "ConjugateTranspose", "NotUnitTriangular", &M, &N, &Alpha, R, &lda, Q, &ldb));
2811:   /* Ainv is Q, overwritten with inverse */

2813: #if defined(PETSC_USE_COMPLEX)
2814:   {
2815:     PetscInt i;
2816:     for (i = 0; i < m * n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
2817:     PetscCall(PetscFree2(A, Ainv));
2818:   }
2819: #endif
2820:   PetscFunctionReturn(PETSC_SUCCESS);
2821: }

2823: /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
2824: static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval, const PetscReal *x, PetscInt ndegree, const PetscInt *degrees, PetscBool Transpose, PetscReal *B)
2825: {
2826:   PetscReal *Bv;
2827:   PetscInt   i, j;

2829:   PetscFunctionBegin;
2830:   PetscCall(PetscMalloc1((ninterval + 1) * ndegree, &Bv));
2831:   /* Point evaluation of L_p on all the source vertices */
2832:   PetscCall(PetscDTLegendreEval(ninterval + 1, x, ndegree, degrees, Bv, NULL, NULL));
2833:   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
2834:   for (i = 0; i < ninterval; i++) {
2835:     for (j = 0; j < ndegree; j++) {
2836:       if (Transpose) B[i + ninterval * j] = Bv[(i + 1) * ndegree + j] - Bv[i * ndegree + j];
2837:       else B[i * ndegree + j] = Bv[(i + 1) * ndegree + j] - Bv[i * ndegree + j];
2838:     }
2839:   }
2840:   PetscCall(PetscFree(Bv));
2841:   PetscFunctionReturn(PETSC_SUCCESS);
2842: }

2844: /*@
2845:   PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals

2847:   Not Collective

2849:   Input Parameters:
2850: + degree  - degree of reconstruction polynomial
2851: . nsource - number of source intervals
2852: . sourcex - sorted coordinates of source cell boundaries (length `nsource`+1)
2853: . ntarget - number of target intervals
2854: - targetx - sorted coordinates of target cell boundaries (length `ntarget`+1)

2856:   Output Parameter:
2857: . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]

2859:   Level: advanced

2861: .seealso: `PetscDTLegendreEval()`
2862: @*/
2863: PetscErrorCode PetscDTReconstructPoly(PetscInt degree, PetscInt nsource, const PetscReal sourcex[], PetscInt ntarget, const PetscReal targetx[], PetscReal R[])
2864: {
2865:   PetscInt     i, j, k, *bdegrees, worksize;
2866:   PetscReal    xmin, xmax, center, hscale, *sourcey, *targety, *Bsource, *Bsinv, *Btarget;
2867:   PetscScalar *tau, *work;

2869:   PetscFunctionBegin;
2870:   PetscAssertPointer(sourcex, 3);
2871:   PetscAssertPointer(targetx, 5);
2872:   PetscAssertPointer(R, 6);
2873:   PetscCheck(degree < nsource, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Reconstruction degree %" PetscInt_FMT " must be less than number of source intervals %" PetscInt_FMT, degree, nsource);
2874:   if (PetscDefined(USE_DEBUG)) {
2875:     for (i = 0; i < nsource; i++) PetscCheck(sourcex[i] < sourcex[i + 1], PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Source interval %" PetscInt_FMT " has negative orientation (%g,%g)", i, (double)sourcex[i], (double)sourcex[i + 1]);
2876:     for (i = 0; i < ntarget; i++) PetscCheck(targetx[i] < targetx[i + 1], PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Target interval %" PetscInt_FMT " has negative orientation (%g,%g)", i, (double)targetx[i], (double)targetx[i + 1]);
2877:   }
2878:   xmin     = PetscMin(sourcex[0], targetx[0]);
2879:   xmax     = PetscMax(sourcex[nsource], targetx[ntarget]);
2880:   center   = (xmin + xmax) / 2;
2881:   hscale   = (xmax - xmin) / 2;
2882:   worksize = nsource;
2883:   PetscCall(PetscMalloc4(degree + 1, &bdegrees, nsource + 1, &sourcey, nsource * (degree + 1), &Bsource, worksize, &work));
2884:   PetscCall(PetscMalloc4(nsource, &tau, nsource * (degree + 1), &Bsinv, ntarget + 1, &targety, ntarget * (degree + 1), &Btarget));
2885:   for (i = 0; i <= nsource; i++) sourcey[i] = (sourcex[i] - center) / hscale;
2886:   for (i = 0; i <= degree; i++) bdegrees[i] = i + 1;
2887:   PetscCall(PetscDTLegendreIntegrate(nsource, sourcey, degree + 1, bdegrees, PETSC_TRUE, Bsource));
2888:   PetscCall(PetscDTPseudoInverseQR(nsource, nsource, degree + 1, Bsource, Bsinv, tau, nsource, work));
2889:   for (i = 0; i <= ntarget; i++) targety[i] = (targetx[i] - center) / hscale;
2890:   PetscCall(PetscDTLegendreIntegrate(ntarget, targety, degree + 1, bdegrees, PETSC_FALSE, Btarget));
2891:   for (i = 0; i < ntarget; i++) {
2892:     PetscReal rowsum = 0;
2893:     for (j = 0; j < nsource; j++) {
2894:       PetscReal sum = 0;
2895:       for (k = 0; k < degree + 1; k++) sum += Btarget[i * (degree + 1) + k] * Bsinv[k * nsource + j];
2896:       R[i * nsource + j] = sum;
2897:       rowsum += sum;
2898:     }
2899:     for (j = 0; j < nsource; j++) R[i * nsource + j] /= rowsum; /* normalize each row */
2900:   }
2901:   PetscCall(PetscFree4(bdegrees, sourcey, Bsource, work));
2902:   PetscCall(PetscFree4(tau, Bsinv, targety, Btarget));
2903:   PetscFunctionReturn(PETSC_SUCCESS);
2904: }

2906: /*@
2907:   PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points

2909:   Not Collective

2911:   Input Parameters:
2912: + n       - the number of GLL nodes
2913: . nodes   - the GLL nodes
2914: . weights - the GLL weights
2915: - f       - the function values at the nodes

2917:   Output Parameter:
2918: . in - the value of the integral

2920:   Level: beginner

2922: .seealso: `PetscDTGaussLobattoLegendreQuadrature()`
2923: @*/
2924: PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n, PetscReal nodes[], PetscReal weights[], const PetscReal f[], PetscReal *in)
2925: {
2926:   PetscInt i;

2928:   PetscFunctionBegin;
2929:   *in = 0.;
2930:   for (i = 0; i < n; i++) *in += f[i] * f[i] * weights[i];
2931:   PetscFunctionReturn(PETSC_SUCCESS);
2932: }

2934: /*@C
2935:   PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element

2937:   Not Collective

2939:   Input Parameters:
2940: + n       - the number of GLL nodes
2941: . nodes   - the GLL nodes, of length `n`
2942: - weights - the GLL weights, of length `n`

2944:   Output Parameter:
2945: . AA - the stiffness element, of size `n` by `n`

2947:   Level: beginner

2949:   Notes:
2950:   Destroy this with `PetscGaussLobattoLegendreElementLaplacianDestroy()`

2952:   You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row-oriented (the array is symmetric)

2954: .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianDestroy()`
2955: @*/
2956: PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n, PetscReal nodes[], PetscReal weights[], PetscReal ***AA)
2957: {
2958:   PetscReal      **A;
2959:   const PetscReal *gllnodes = nodes;
2960:   const PetscInt   p        = n - 1;
2961:   PetscReal        z0, z1, z2 = -1, x, Lpj, Lpr;
2962:   PetscInt         i, j, nn, r;

2964:   PetscFunctionBegin;
2965:   PetscCall(PetscMalloc1(n, &A));
2966:   PetscCall(PetscMalloc1(n * n, &A[0]));
2967:   for (i = 1; i < n; i++) A[i] = A[i - 1] + n;

2969:   for (j = 1; j < p; j++) {
2970:     x  = gllnodes[j];
2971:     z0 = 1.;
2972:     z1 = x;
2973:     for (nn = 1; nn < p; nn++) {
2974:       z2 = x * z1 * (2. * ((PetscReal)nn) + 1.) / (((PetscReal)nn) + 1.) - z0 * (((PetscReal)nn) / (((PetscReal)nn) + 1.));
2975:       z0 = z1;
2976:       z1 = z2;
2977:     }
2978:     Lpj = z2;
2979:     for (r = 1; r < p; r++) {
2980:       if (r == j) {
2981:         A[j][j] = 2. / (3. * (1. - gllnodes[j] * gllnodes[j]) * Lpj * Lpj);
2982:       } else {
2983:         x  = gllnodes[r];
2984:         z0 = 1.;
2985:         z1 = x;
2986:         for (nn = 1; nn < p; nn++) {
2987:           z2 = x * z1 * (2. * ((PetscReal)nn) + 1.) / (((PetscReal)nn) + 1.) - z0 * (((PetscReal)nn) / (((PetscReal)nn) + 1.));
2988:           z0 = z1;
2989:           z1 = z2;
2990:         }
2991:         Lpr     = z2;
2992:         A[r][j] = 4. / (((PetscReal)p) * (((PetscReal)p) + 1.) * Lpj * Lpr * (gllnodes[j] - gllnodes[r]) * (gllnodes[j] - gllnodes[r]));
2993:       }
2994:     }
2995:   }
2996:   for (j = 1; j < p + 1; j++) {
2997:     x  = gllnodes[j];
2998:     z0 = 1.;
2999:     z1 = x;
3000:     for (nn = 1; nn < p; nn++) {
3001:       z2 = x * z1 * (2. * ((PetscReal)nn) + 1.) / (((PetscReal)nn) + 1.) - z0 * (((PetscReal)nn) / (((PetscReal)nn) + 1.));
3002:       z0 = z1;
3003:       z1 = z2;
3004:     }
3005:     Lpj     = z2;
3006:     A[j][0] = 4. * PetscPowRealInt(-1., p) / (((PetscReal)p) * (((PetscReal)p) + 1.) * Lpj * (1. + gllnodes[j]) * (1. + gllnodes[j]));
3007:     A[0][j] = A[j][0];
3008:   }
3009:   for (j = 0; j < p; j++) {
3010:     x  = gllnodes[j];
3011:     z0 = 1.;
3012:     z1 = x;
3013:     for (nn = 1; nn < p; nn++) {
3014:       z2 = x * z1 * (2. * ((PetscReal)nn) + 1.) / (((PetscReal)nn) + 1.) - z0 * (((PetscReal)nn) / (((PetscReal)nn) + 1.));
3015:       z0 = z1;
3016:       z1 = z2;
3017:     }
3018:     Lpj = z2;

3020:     A[p][j] = 4. / (((PetscReal)p) * (((PetscReal)p) + 1.) * Lpj * (1. - gllnodes[j]) * (1. - gllnodes[j]));
3021:     A[j][p] = A[p][j];
3022:   }
3023:   A[0][0] = 0.5 + (((PetscReal)p) * (((PetscReal)p) + 1.) - 2.) / 6.;
3024:   A[p][p] = A[0][0];
3025:   *AA     = A;
3026:   PetscFunctionReturn(PETSC_SUCCESS);
3027: }

3029: /*@C
3030:   PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element created with `PetscGaussLobattoLegendreElementLaplacianCreate()`

3032:   Not Collective

3034:   Input Parameters:
3035: + n       - the number of GLL nodes
3036: . nodes   - the GLL nodes, ignored
3037: . weights - the GLL weightss, ignored
3038: - AA      - the stiffness element from `PetscGaussLobattoLegendreElementLaplacianCreate()`

3040:   Level: beginner

3042: .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianCreate()`
3043: @*/
3044: PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n, PetscReal nodes[], PetscReal weights[], PetscReal ***AA)
3045: {
3046:   PetscFunctionBegin;
3047:   PetscCall(PetscFree((*AA)[0]));
3048:   PetscCall(PetscFree(*AA));
3049:   *AA = NULL;
3050:   PetscFunctionReturn(PETSC_SUCCESS);
3051: }

3053: /*@C
3054:   PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element

3056:   Not Collective

3058:   Input Parameters:
3059: + n       - the number of GLL nodes
3060: . nodes   - the GLL nodes, of length `n`
3061: - weights - the GLL weights, of length `n`

3063:   Output Parameters:
3064: + AA  - the stiffness element, of dimension `n` by `n`
3065: - AAT - the transpose of AA (pass in `NULL` if you do not need this array), of dimension `n` by `n`

3067:   Level: beginner

3069:   Notes:
3070:   Destroy this with `PetscGaussLobattoLegendreElementGradientDestroy()`

3072:   You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row-oriented

3074: .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianDestroy()`, `PetscGaussLobattoLegendreElementGradientDestroy()`
3075: @*/
3076: PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n, PetscReal nodes[], PetscReal weights[], PetscReal ***AA, PetscReal ***AAT)
3077: {
3078:   PetscReal      **A, **AT = NULL;
3079:   const PetscReal *gllnodes = nodes;
3080:   const PetscInt   p        = n - 1;
3081:   PetscReal        Li, Lj, d0;
3082:   PetscInt         i, j;

3084:   PetscFunctionBegin;
3085:   PetscCall(PetscMalloc1(n, &A));
3086:   PetscCall(PetscMalloc1(n * n, &A[0]));
3087:   for (i = 1; i < n; i++) A[i] = A[i - 1] + n;

3089:   if (AAT) {
3090:     PetscCall(PetscMalloc1(n, &AT));
3091:     PetscCall(PetscMalloc1(n * n, &AT[0]));
3092:     for (i = 1; i < n; i++) AT[i] = AT[i - 1] + n;
3093:   }

3095:   if (n == 1) A[0][0] = 0.;
3096:   d0 = (PetscReal)p * ((PetscReal)p + 1.) / 4.;
3097:   for (i = 0; i < n; i++) {
3098:     for (j = 0; j < n; j++) {
3099:       A[i][j] = 0.;
3100:       PetscCall(PetscDTComputeJacobi(0., 0., p, gllnodes[i], &Li));
3101:       PetscCall(PetscDTComputeJacobi(0., 0., p, gllnodes[j], &Lj));
3102:       if (i != j) A[i][j] = Li / (Lj * (gllnodes[i] - gllnodes[j]));
3103:       if ((j == i) && (i == 0)) A[i][j] = -d0;
3104:       if (j == i && i == p) A[i][j] = d0;
3105:       if (AT) AT[j][i] = A[i][j];
3106:     }
3107:   }
3108:   if (AAT) *AAT = AT;
3109:   *AA = A;
3110:   PetscFunctionReturn(PETSC_SUCCESS);
3111: }

3113: /*@C
3114:   PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with `PetscGaussLobattoLegendreElementGradientCreate()`

3116:   Not Collective

3118:   Input Parameters:
3119: + n       - the number of GLL nodes
3120: . nodes   - the GLL nodes, ignored
3121: . weights - the GLL weights, ignored
3122: . AA      - the stiffness element obtained with `PetscGaussLobattoLegendreElementGradientCreate()`
3123: - AAT     - the transpose of the element obtained with `PetscGaussLobattoLegendreElementGradientCreate()`

3125:   Level: beginner

3127: .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianCreate()`, `PetscGaussLobattoLegendreElementAdvectionCreate()`
3128: @*/
3129: PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n, PetscReal nodes[], PetscReal weights[], PetscReal ***AA, PetscReal ***AAT)
3130: {
3131:   PetscFunctionBegin;
3132:   PetscCall(PetscFree((*AA)[0]));
3133:   PetscCall(PetscFree(*AA));
3134:   *AA = NULL;
3135:   if (AAT) {
3136:     PetscCall(PetscFree((*AAT)[0]));
3137:     PetscCall(PetscFree(*AAT));
3138:     *AAT = NULL;
3139:   }
3140:   PetscFunctionReturn(PETSC_SUCCESS);
3141: }

3143: /*@C
3144:   PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element

3146:   Not Collective

3148:   Input Parameters:
3149: + n       - the number of GLL nodes
3150: . nodes   - the GLL nodes, of length `n`
3151: - weights - the GLL weights, of length `n`

3153:   Output Parameter:
3154: . AA - the stiffness element, of dimension `n` by `n`

3156:   Level: beginner

3158:   Notes:
3159:   Destroy this with `PetscGaussLobattoLegendreElementAdvectionDestroy()`

3161:   This is the same as the Gradient operator multiplied by the diagonal mass matrix

3163:   You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row-oriented

3165: .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianCreate()`, `PetscGaussLobattoLegendreElementAdvectionDestroy()`
3166: @*/
3167: PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n, PetscReal nodes[], PetscReal weights[], PetscReal ***AA)
3168: {
3169:   PetscReal      **D;
3170:   const PetscReal *gllweights = weights;
3171:   const PetscInt   glln       = n;
3172:   PetscInt         j;

3174:   PetscFunctionBegin;
3175:   PetscCall(PetscGaussLobattoLegendreElementGradientCreate(n, nodes, weights, &D, NULL));
3176:   for (PetscInt i = 0; i < glln; i++) {
3177:     for (j = 0; j < glln; j++) D[i][j] = gllweights[i] * D[i][j];
3178:   }
3179:   *AA = D;
3180:   PetscFunctionReturn(PETSC_SUCCESS);
3181: }

3183: /*@C
3184:   PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element created with `PetscGaussLobattoLegendreElementAdvectionCreate()`

3186:   Not Collective

3188:   Input Parameters:
3189: + n       - the number of GLL nodes
3190: . nodes   - the GLL nodes, ignored
3191: . weights - the GLL weights, ignored
3192: - AA      - advection obtained with `PetscGaussLobattoLegendreElementAdvectionCreate()`

3194:   Level: beginner

3196: .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementAdvectionCreate()`
3197: @*/
3198: PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n, PetscReal nodes[], PetscReal weights[], PetscReal ***AA)
3199: {
3200:   PetscFunctionBegin;
3201:   PetscCall(PetscFree((*AA)[0]));
3202:   PetscCall(PetscFree(*AA));
3203:   *AA = NULL;
3204:   PetscFunctionReturn(PETSC_SUCCESS);
3205: }

3207: PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA)
3208: {
3209:   PetscReal      **A;
3210:   const PetscReal *gllweights = weights;
3211:   const PetscInt   glln       = n;

3213:   PetscFunctionBegin;
3214:   PetscCall(PetscMalloc1(glln, &A));
3215:   PetscCall(PetscMalloc1(glln * glln, &A[0]));
3216:   for (PetscInt i = 1; i < glln; i++) A[i] = A[i - 1] + glln;
3217:   if (glln == 1) A[0][0] = 0.;
3218:   for (PetscInt i = 0; i < glln; i++) {
3219:     for (PetscInt j = 0; j < glln; j++) {
3220:       A[i][j] = 0.;
3221:       if (j == i) A[i][j] = gllweights[i];
3222:     }
3223:   }
3224:   *AA = A;
3225:   PetscFunctionReturn(PETSC_SUCCESS);
3226: }

3228: PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA)
3229: {
3230:   PetscFunctionBegin;
3231:   PetscCall(PetscFree((*AA)[0]));
3232:   PetscCall(PetscFree(*AA));
3233:   *AA = NULL;
3234:   PetscFunctionReturn(PETSC_SUCCESS);
3235: }

3237: /*@
3238:   PetscDTIndexToBary - convert an index into a barycentric coordinate.

3240:   Input Parameters:
3241: + len   - the desired length of the barycentric tuple (usually 1 more than the dimension it represents, so a barycentric coordinate in a triangle has length 3)
3242: . sum   - the value that the sum of the barycentric coordinates (which will be non-negative integers) should sum to
3243: - index - the index to convert: should be >= 0 and < Binomial(len - 1 + sum, sum)

3245:   Output Parameter:
3246: . coord - will be filled with the barycentric coordinate, of length `n`

3248:   Level: beginner

3250:   Note:
3251:   The indices map to barycentric coordinates in lexicographic order, where the first index is the
3252:   least significant and the last index is the most significant.

3254: .seealso: `PetscDTBaryToIndex()`
3255: @*/
3256: PetscErrorCode PetscDTIndexToBary(PetscInt len, PetscInt sum, PetscInt index, PetscInt coord[])
3257: {
3258:   PetscInt c, d, s, total, subtotal, nexttotal;

3260:   PetscFunctionBeginHot;
3261:   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
3262:   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index must be non-negative");
3263:   if (!len) {
3264:     PetscCheck(!sum && !index, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid index or sum for length 0 barycentric coordinate");
3265:     PetscFunctionReturn(PETSC_SUCCESS);
3266:   }
3267:   for (c = 1, total = 1; c <= len; c++) {
3268:     /* total is the number of ways to have a tuple of length c with sum */
3269:     if (index < total) break;
3270:     total = (total * (sum + c)) / c;
3271:   }
3272:   PetscCheck(c <= len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index out of range");
3273:   for (d = c; d < len; d++) coord[d] = 0;
3274:   for (s = 0, subtotal = 1, nexttotal = 1; c > 0;) {
3275:     /* subtotal is the number of ways to have a tuple of length c with sum s */
3276:     /* nexttotal is the number of ways to have a tuple of length c-1 with sum s */
3277:     if ((index + subtotal) >= total) {
3278:       coord[--c] = sum - s;
3279:       index -= (total - subtotal);
3280:       sum       = s;
3281:       total     = nexttotal;
3282:       subtotal  = 1;
3283:       nexttotal = 1;
3284:       s         = 0;
3285:     } else {
3286:       subtotal  = (subtotal * (c + s)) / (s + 1);
3287:       nexttotal = (nexttotal * (c - 1 + s)) / (s + 1);
3288:       s++;
3289:     }
3290:   }
3291:   PetscFunctionReturn(PETSC_SUCCESS);
3292: }

3294: /*@
3295:   PetscDTBaryToIndex - convert a barycentric coordinate to an index

3297:   Input Parameters:
3298: + len   - the desired length of the barycentric tuple (usually 1 more than the dimension it represents, so a barycentric coordinate in a triangle has length 3)
3299: . sum   - the value that the sum of the barycentric coordinates (which will be non-negative integers) should sum to
3300: - coord - a barycentric coordinate with the given length `len` and `sum`

3302:   Output Parameter:
3303: . index - the unique index for the coordinate, >= 0 and < Binomial(len - 1 + sum, sum)

3305:   Level: beginner

3307:   Note:
3308:   The indices map to barycentric coordinates in lexicographic order, where the first index is the
3309:   least significant and the last index is the most significant.

3311: .seealso: `PetscDTIndexToBary`
3312: @*/
3313: PetscErrorCode PetscDTBaryToIndex(PetscInt len, PetscInt sum, const PetscInt coord[], PetscInt *index)
3314: {
3315:   PetscInt c;
3316:   PetscInt i;
3317:   PetscInt total;

3319:   PetscFunctionBeginHot;
3320:   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
3321:   if (!len) {
3322:     PetscCheck(!sum, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid index or sum for length 0 barycentric coordinate");
3323:     *index = 0;
3324:     PetscFunctionReturn(PETSC_SUCCESS);
3325:   }
3326:   for (c = 1, total = 1; c < len; c++) total = (total * (sum + c)) / c;
3327:   i = total - 1;
3328:   c = len - 1;
3329:   sum -= coord[c];
3330:   while (sum > 0) {
3331:     PetscInt subtotal;
3332:     PetscInt s;

3334:     for (s = 1, subtotal = 1; s < sum; s++) subtotal = (subtotal * (c + s)) / s;
3335:     i -= subtotal;
3336:     sum -= coord[--c];
3337:   }
3338:   *index = i;
3339:   PetscFunctionReturn(PETSC_SUCCESS);
3340: }

3342: /*@
3343:   PetscQuadratureComputePermutations - Compute permutations of quadrature points corresponding to domain orientations

3345:   Input Parameter:
3346: . quad - The `PetscQuadrature`

3348:   Output Parameters:
3349: + Np   - The number of domain orientations
3350: - perm - An array of `IS` permutations, one for ech orientation,

3352:   Level: developer

3354: .seealso: `PetscQuadratureSetCellType()`, `PetscQuadrature`
3355: @*/
3356: PetscErrorCode PetscQuadratureComputePermutations(PetscQuadrature quad, PeOp PetscInt *Np, IS *perm[])
3357: {
3358:   DMPolytopeType   ct;
3359:   const PetscReal *xq, *wq;
3360:   PetscInt         dim, qdim, d, Na, o, Nq, q, qp;

3362:   PetscFunctionBegin;
3363:   PetscCall(PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &xq, &wq));
3364:   PetscCall(PetscQuadratureGetCellType(quad, &ct));
3365:   dim = DMPolytopeTypeGetDim(ct);
3366:   Na  = DMPolytopeTypeGetNumArrangements(ct);
3367:   PetscCall(PetscMalloc1(Na, perm));
3368:   if (Np) *Np = Na;
3369:   Na /= 2;
3370:   for (o = -Na; o < Na; ++o) {
3371:     DM        refdm;
3372:     PetscInt *idx;
3373:     PetscReal xi0[3] = {-1., -1., -1.}, v0[3], J[9], detJ, txq[3];
3374:     PetscBool flg;

3376:     PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &refdm));
3377:     PetscCall(DMPlexOrientPoint(refdm, 0, o));
3378:     PetscCall(DMPlexComputeCellGeometryFEM(refdm, 0, NULL, v0, J, NULL, &detJ));
3379:     PetscCall(PetscMalloc1(Nq, &idx));
3380:     for (q = 0; q < Nq; ++q) {
3381:       CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], txq);
3382:       for (qp = 0; qp < Nq; ++qp) {
3383:         PetscReal diff = 0.;

3385:         for (d = 0; d < dim; ++d) diff += PetscAbsReal(txq[d] - xq[qp * dim + d]);
3386:         if (diff < PETSC_SMALL) break;
3387:       }
3388:       PetscCheck(qp < Nq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Transformed quad point %" PetscInt_FMT " does not match another quad point", q);
3389:       idx[q] = qp;
3390:     }
3391:     PetscCall(DMDestroy(&refdm));
3392:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nq, idx, PETSC_OWN_POINTER, &(*perm)[o + Na]));
3393:     PetscCall(ISGetInfo((*perm)[o + Na], IS_PERMUTATION, IS_LOCAL, PETSC_TRUE, &flg));
3394:     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ordering for orientation %" PetscInt_FMT " was not a permutation", o);
3395:     PetscCall(ISSetPermutation((*perm)[o + Na]));
3396:   }
3397:   if (!Na) (*perm)[0] = NULL;
3398:   PetscFunctionReturn(PETSC_SUCCESS);
3399: }

3401: /*@
3402:   PetscDTCreateQuadratureByCell - Create default quadrature for a given cell

3404:   Not collective

3406:   Input Parameters:
3407: + ct     - The integration domain
3408: . qorder - The desired quadrature order
3409: - qtype  - The type of simplex quadrature, or PETSCDTSIMPLEXQUAD_DEFAULT

3411:   Output Parameters:
3412: + q  - The cell quadrature
3413: - fq - The face quadrature

3415:   Level: developer

3417: .seealso: `PetscDTCreateDefaultQuadrature()`, `PetscFECreateDefault()`, `PetscDTGaussTensorQuadrature()`, `PetscDTSimplexQuadrature()`, `PetscDTTensorQuadratureCreate()`
3418: @*/
3419: PetscErrorCode PetscDTCreateQuadratureByCell(DMPolytopeType ct, PetscInt qorder, PetscDTSimplexQuadratureType qtype, PetscQuadrature *q, PetscQuadrature *fq)
3420: {
3421:   const PetscInt quadPointsPerEdge = PetscMax(qorder + 1, 1);
3422:   const PetscInt dim               = DMPolytopeTypeGetDim(ct);

3424:   PetscFunctionBegin;
3425:   switch (ct) {
3426:   case DM_POLYTOPE_POINT: {
3427:     PetscReal *x, *w;

3429:     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
3430:     PetscCall(PetscQuadratureSetCellType(*q, ct));
3431:     PetscCall(PetscQuadratureSetOrder(*q, 0));
3432:     PetscCall(PetscMalloc1(1, &x));
3433:     PetscCall(PetscMalloc1(1, &w));
3434:     x[0] = 0.;
3435:     w[0] = 1.;
3436:     PetscCall(PetscQuadratureSetData(*q, dim, 1, 1, x, w));
3437:     *fq = NULL;
3438:   } break;
3439:   case DM_POLYTOPE_SEGMENT:
3440:   case DM_POLYTOPE_POINT_PRISM_TENSOR:
3441:   case DM_POLYTOPE_QUADRILATERAL:
3442:   case DM_POLYTOPE_SEG_PRISM_TENSOR:
3443:   case DM_POLYTOPE_HEXAHEDRON:
3444:   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
3445:     PetscCall(PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, q));
3446:     PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, fq));
3447:     break;
3448:   case DM_POLYTOPE_TRIANGLE:
3449:   case DM_POLYTOPE_TETRAHEDRON:
3450:     PetscCall(PetscDTSimplexQuadrature(dim, 2 * qorder, qtype, q));
3451:     PetscCall(PetscDTSimplexQuadrature(dim - 1, 2 * qorder, qtype, fq));
3452:     break;
3453:   case DM_POLYTOPE_TRI_PRISM:
3454:   case DM_POLYTOPE_TRI_PRISM_TENSOR: {
3455:     PetscQuadrature q1, q2;

3457:     // TODO: this should be able to use symmetric rules, but doing so causes tests to fail
3458:     PetscCall(PetscDTSimplexQuadrature(2, 2 * qorder, PETSCDTSIMPLEXQUAD_CONIC, &q1));
3459:     PetscCall(PetscDTGaussTensorQuadrature(1, 1, quadPointsPerEdge, -1.0, 1.0, &q2));
3460:     PetscCall(PetscDTTensorQuadratureCreate(q1, q2, q));
3461:     PetscCall(PetscQuadratureDestroy(&q2));
3462:     *fq = q1;
3463:     /* TODO Need separate quadratures for each face */
3464:   } break;
3465:   default:
3466:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No quadrature for celltype %s", DMPolytopeTypes[PetscMin(ct, DM_POLYTOPE_UNKNOWN)]);
3467:   }
3468:   PetscFunctionReturn(PETSC_SUCCESS);
3469: }

3471: /*@
3472:   PetscDTCreateDefaultQuadrature - Create default quadrature for a given cell

3474:   Not collective

3476:   Input Parameters:
3477: + ct     - The integration domain
3478: - qorder - The desired quadrature order

3480:   Output Parameters:
3481: + q  - The cell quadrature
3482: - fq - The face quadrature

3484:   Level: developer

3486: .seealso: `PetscDTCreateQuadratureByCell()`, `PetscFECreateDefault()`, `PetscDTGaussTensorQuadrature()`, `PetscDTSimplexQuadrature()`, `PetscDTTensorQuadratureCreate()`
3487: @*/
3488: PetscErrorCode PetscDTCreateDefaultQuadrature(DMPolytopeType ct, PetscInt qorder, PetscQuadrature *q, PetscQuadrature *fq)
3489: {
3490:   PetscFunctionBegin;
3491:   PetscCall(PetscDTCreateQuadratureByCell(ct, qorder, PETSCDTSIMPLEXQUAD_DEFAULT, q, fq));
3492:   PetscFunctionReturn(PETSC_SUCCESS);
3493: }