xref: /petsc/src/dm/dt/interface/dt.c (revision 07218a29b4d399ea5732e51342b017bba1a66494)
137045ce4SJed Brown /* Discretization tools */
237045ce4SJed Brown 
30c35b76eSJed Brown #include <petscdt.h> /*I "petscdt.h" I*/
437045ce4SJed Brown #include <petscblaslapack.h>
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
6af0996ceSBarry Smith #include <petsc/private/dtimpl.h>
7*07218a29SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */
8665c2dedSJed Brown #include <petscviewer.h>
959804f93SMatthew G. Knepley #include <petscdmplex.h>
1059804f93SMatthew G. Knepley #include <petscdmshell.h>
1137045ce4SJed Brown 
1298c04793SMatthew G. Knepley #if defined(PETSC_HAVE_MPFR)
1398c04793SMatthew G. Knepley   #include <mpfr.h>
1498c04793SMatthew G. Knepley #endif
1598c04793SMatthew G. Knepley 
16d3c69ad0SToby Isaac const char *const        PetscDTNodeTypes_shifted[] = {"default", "gaussjacobi", "equispaced", "tanhsinh", "PETSCDTNODES_", NULL};
17d3c69ad0SToby Isaac const char *const *const PetscDTNodeTypes           = PetscDTNodeTypes_shifted + 1;
18d3c69ad0SToby Isaac 
19d3c69ad0SToby Isaac const char *const        PetscDTSimplexQuadratureTypes_shifted[] = {"default", "conic", "minsym", "PETSCDTSIMPLEXQUAD_", NULL};
20d3c69ad0SToby Isaac const char *const *const PetscDTSimplexQuadratureTypes           = PetscDTSimplexQuadratureTypes_shifted + 1;
21d4afb720SToby Isaac 
22e6a796c3SToby Isaac static PetscBool GolubWelschCite       = PETSC_FALSE;
23e6a796c3SToby Isaac const char       GolubWelschCitation[] = "@article{GolubWelsch1969,\n"
240bfcf5a5SMatthew G. Knepley                                          "  author  = {Golub and Welsch},\n"
250bfcf5a5SMatthew G. Knepley                                          "  title   = {Calculation of Quadrature Rules},\n"
260bfcf5a5SMatthew G. Knepley                                          "  journal = {Math. Comp.},\n"
270bfcf5a5SMatthew G. Knepley                                          "  volume  = {23},\n"
280bfcf5a5SMatthew G. Knepley                                          "  number  = {106},\n"
290bfcf5a5SMatthew G. Knepley                                          "  pages   = {221--230},\n"
300bfcf5a5SMatthew G. Knepley                                          "  year    = {1969}\n}\n";
310bfcf5a5SMatthew G. Knepley 
32c4762a1bSJed Brown /* Numerical tests in src/dm/dt/tests/ex1.c show that when computing the nodes and weights of Gauss-Jacobi
3394e21283SToby Isaac    quadrature rules:
34e6a796c3SToby Isaac 
3594e21283SToby Isaac    - in double precision, Newton's method and Golub & Welsch both work for moderate degrees (< 100),
3694e21283SToby Isaac    - in single precision, Newton's method starts producing incorrect roots around n = 15, but
3794e21283SToby Isaac      the weights from Golub & Welsch become a problem before then: they produces errors
3894e21283SToby Isaac      in computing the Jacobi-polynomial Gram matrix around n = 6.
3994e21283SToby Isaac 
4094e21283SToby Isaac    So we default to Newton's method (required fewer dependencies) */
4194e21283SToby Isaac PetscBool PetscDTGaussQuadratureNewton_Internal = PETSC_TRUE;
422cd22861SMatthew G. Knepley 
432cd22861SMatthew G. Knepley PetscClassId PETSCQUADRATURE_CLASSID = 0;
442cd22861SMatthew G. Knepley 
4540d8ff71SMatthew G. Knepley /*@
46dce8aebaSBarry Smith   PetscQuadratureCreate - Create a `PetscQuadrature` object
4740d8ff71SMatthew G. Knepley 
48d083f849SBarry Smith   Collective
4940d8ff71SMatthew G. Knepley 
5040d8ff71SMatthew G. Knepley   Input Parameter:
51dce8aebaSBarry Smith . comm - The communicator for the `PetscQuadrature` object
5240d8ff71SMatthew G. Knepley 
5340d8ff71SMatthew G. Knepley   Output Parameter:
5420f4b53cSBarry Smith . q  - The `PetscQuadrature` object
5540d8ff71SMatthew G. Knepley 
5640d8ff71SMatthew G. Knepley   Level: beginner
5740d8ff71SMatthew G. Knepley 
58dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `Petscquadraturedestroy()`, `PetscQuadratureGetData()`
5940d8ff71SMatthew G. Knepley @*/
60d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
61d71ae5a4SJacob Faibussowitsch {
6221454ff5SMatthew G. Knepley   PetscFunctionBegin;
6321454ff5SMatthew G. Knepley   PetscValidPointer(q, 2);
649566063dSJacob Faibussowitsch   PetscCall(DMInitializePackage());
659566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*q, PETSCQUADRATURE_CLASSID, "PetscQuadrature", "Quadrature", "DT", comm, PetscQuadratureDestroy, PetscQuadratureView));
6621454ff5SMatthew G. Knepley   (*q)->dim       = -1;
67a6b92713SMatthew G. Knepley   (*q)->Nc        = 1;
68bcede257SMatthew G. Knepley   (*q)->order     = -1;
6921454ff5SMatthew G. Knepley   (*q)->numPoints = 0;
7021454ff5SMatthew G. Knepley   (*q)->points    = NULL;
7121454ff5SMatthew G. Knepley   (*q)->weights   = NULL;
723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7321454ff5SMatthew G. Knepley }
7421454ff5SMatthew G. Knepley 
75c9638911SMatthew G. Knepley /*@
76dce8aebaSBarry Smith   PetscQuadratureDuplicate - Create a deep copy of the `PetscQuadrature` object
77c9638911SMatthew G. Knepley 
7820f4b53cSBarry Smith   Collective
79c9638911SMatthew G. Knepley 
80c9638911SMatthew G. Knepley   Input Parameter:
81dce8aebaSBarry Smith . q  - The `PetscQuadrature` object
82c9638911SMatthew G. Knepley 
83c9638911SMatthew G. Knepley   Output Parameter:
84dce8aebaSBarry Smith . r  - The new `PetscQuadrature` object
85c9638911SMatthew G. Knepley 
86c9638911SMatthew G. Knepley   Level: beginner
87c9638911SMatthew G. Knepley 
88dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`, `PetscQuadratureDestroy()`, `PetscQuadratureGetData()`
89c9638911SMatthew G. Knepley @*/
90d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
91d71ae5a4SJacob Faibussowitsch {
92a6b92713SMatthew G. Knepley   PetscInt         order, dim, Nc, Nq;
93c9638911SMatthew G. Knepley   const PetscReal *points, *weights;
94c9638911SMatthew G. Knepley   PetscReal       *p, *w;
95c9638911SMatthew G. Knepley 
96c9638911SMatthew G. Knepley   PetscFunctionBegin;
97064a246eSJacob Faibussowitsch   PetscValidPointer(q, 1);
989566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)q), r));
999566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetOrder(q, &order));
1009566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetOrder(*r, order));
1019566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights));
1029566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nq * dim, &p));
1039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nq * Nc, &w));
1049566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(p, points, Nq * dim));
1059566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(w, weights, Nc * Nq));
1069566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*r, dim, Nc, Nq, p, w));
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108c9638911SMatthew G. Knepley }
109c9638911SMatthew G. Knepley 
11040d8ff71SMatthew G. Knepley /*@
111dce8aebaSBarry Smith   PetscQuadratureDestroy - Destroys a `PetscQuadrature` object
11240d8ff71SMatthew G. Knepley 
11320f4b53cSBarry Smith   Collective
11440d8ff71SMatthew G. Knepley 
11540d8ff71SMatthew G. Knepley   Input Parameter:
116dce8aebaSBarry Smith . q  - The `PetscQuadrature` object
11740d8ff71SMatthew G. Knepley 
11840d8ff71SMatthew G. Knepley   Level: beginner
11940d8ff71SMatthew G. Knepley 
120dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`, `PetscQuadratureGetData()`
12140d8ff71SMatthew G. Knepley @*/
122d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
123d71ae5a4SJacob Faibussowitsch {
124bfa639d9SMatthew G. Knepley   PetscFunctionBegin;
1253ba16761SJacob Faibussowitsch   if (!*q) PetscFunctionReturn(PETSC_SUCCESS);
1262cd22861SMatthew G. Knepley   PetscValidHeaderSpecific((*q), PETSCQUADRATURE_CLASSID, 1);
12721454ff5SMatthew G. Knepley   if (--((PetscObject)(*q))->refct > 0) {
12821454ff5SMatthew G. Knepley     *q = NULL;
1293ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
13021454ff5SMatthew G. Knepley   }
1319566063dSJacob Faibussowitsch   PetscCall(PetscFree((*q)->points));
1329566063dSJacob Faibussowitsch   PetscCall(PetscFree((*q)->weights));
1339566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(q));
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13521454ff5SMatthew G. Knepley }
13621454ff5SMatthew G. Knepley 
137bcede257SMatthew G. Knepley /*@
138dce8aebaSBarry Smith   PetscQuadratureGetOrder - Return the order of the method in the `PetscQuadrature`
139bcede257SMatthew G. Knepley 
14020f4b53cSBarry Smith   Not Collective
141bcede257SMatthew G. Knepley 
142bcede257SMatthew G. Knepley   Input Parameter:
143dce8aebaSBarry Smith . q - The `PetscQuadrature` object
144bcede257SMatthew G. Knepley 
145bcede257SMatthew G. Knepley   Output Parameter:
146bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
147bcede257SMatthew G. Knepley 
148bcede257SMatthew G. Knepley   Level: intermediate
149bcede257SMatthew G. Knepley 
150dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureSetOrder()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
151bcede257SMatthew G. Knepley @*/
152d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
153d71ae5a4SJacob Faibussowitsch {
154bcede257SMatthew G. Knepley   PetscFunctionBegin;
1552cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
156dadcf809SJacob Faibussowitsch   PetscValidIntPointer(order, 2);
157bcede257SMatthew G. Knepley   *order = q->order;
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159bcede257SMatthew G. Knepley }
160bcede257SMatthew G. Knepley 
161bcede257SMatthew G. Knepley /*@
162dce8aebaSBarry Smith   PetscQuadratureSetOrder - Set the order of the method in the `PetscQuadrature`
163bcede257SMatthew G. Knepley 
16420f4b53cSBarry Smith   Not Collective
165bcede257SMatthew G. Knepley 
166bcede257SMatthew G. Knepley   Input Parameters:
167dce8aebaSBarry Smith + q - The `PetscQuadrature` object
168bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
169bcede257SMatthew G. Knepley 
170bcede257SMatthew G. Knepley   Level: intermediate
171bcede257SMatthew G. Knepley 
172dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureGetOrder()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
173bcede257SMatthew G. Knepley @*/
174d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
175d71ae5a4SJacob Faibussowitsch {
176bcede257SMatthew G. Knepley   PetscFunctionBegin;
1772cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
178bcede257SMatthew G. Knepley   q->order = order;
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
180bcede257SMatthew G. Knepley }
181bcede257SMatthew G. Knepley 
182a6b92713SMatthew G. Knepley /*@
183a6b92713SMatthew G. Knepley   PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated
184a6b92713SMatthew G. Knepley 
18520f4b53cSBarry Smith   Not Collective
186a6b92713SMatthew G. Knepley 
187a6b92713SMatthew G. Knepley   Input Parameter:
188dce8aebaSBarry Smith . q - The `PetscQuadrature` object
189a6b92713SMatthew G. Knepley 
190a6b92713SMatthew G. Knepley   Output Parameter:
191a6b92713SMatthew G. Knepley . Nc - The number of components
192a6b92713SMatthew G. Knepley 
19320f4b53cSBarry Smith   Level: intermediate
19420f4b53cSBarry Smith 
195dce8aebaSBarry Smith   Note:
196dce8aebaSBarry Smith   We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
197a6b92713SMatthew G. Knepley 
198dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureSetNumComponents()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
199a6b92713SMatthew G. Knepley @*/
200d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc)
201d71ae5a4SJacob Faibussowitsch {
202a6b92713SMatthew G. Knepley   PetscFunctionBegin;
2032cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
204dadcf809SJacob Faibussowitsch   PetscValidIntPointer(Nc, 2);
205a6b92713SMatthew G. Knepley   *Nc = q->Nc;
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207a6b92713SMatthew G. Knepley }
208a6b92713SMatthew G. Knepley 
209a6b92713SMatthew G. Knepley /*@
210a6b92713SMatthew G. Knepley   PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated
211a6b92713SMatthew G. Knepley 
21220f4b53cSBarry Smith   Not Collective
213a6b92713SMatthew G. Knepley 
214a6b92713SMatthew G. Knepley   Input Parameters:
215a6b92713SMatthew G. Knepley + q  - The PetscQuadrature object
216a6b92713SMatthew G. Knepley - Nc - The number of components
217a6b92713SMatthew G. Knepley 
21820f4b53cSBarry Smith   Level: intermediate
21920f4b53cSBarry Smith 
220dce8aebaSBarry Smith   Note:
221dce8aebaSBarry Smith   We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
222a6b92713SMatthew G. Knepley 
223dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureGetNumComponents()`, `PetscQuadratureGetData()`, `PetscQuadratureSetData()`
224a6b92713SMatthew G. Knepley @*/
225d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc)
226d71ae5a4SJacob Faibussowitsch {
227a6b92713SMatthew G. Knepley   PetscFunctionBegin;
2282cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
229a6b92713SMatthew G. Knepley   q->Nc = Nc;
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
231a6b92713SMatthew G. Knepley }
232a6b92713SMatthew G. Knepley 
23340d8ff71SMatthew G. Knepley /*@C
234dce8aebaSBarry Smith   PetscQuadratureGetData - Returns the data defining the `PetscQuadrature`
23540d8ff71SMatthew G. Knepley 
23620f4b53cSBarry Smith   Not Collective
23740d8ff71SMatthew G. Knepley 
23840d8ff71SMatthew G. Knepley   Input Parameter:
239dce8aebaSBarry Smith . q  - The `PetscQuadrature` object
24040d8ff71SMatthew G. Knepley 
24140d8ff71SMatthew G. Knepley   Output Parameters:
24240d8ff71SMatthew G. Knepley + dim - The spatial dimension
243805e7170SToby Isaac . Nc - The number of components
24440d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
24540d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
24640d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
24740d8ff71SMatthew G. Knepley 
24840d8ff71SMatthew G. Knepley   Level: intermediate
24940d8ff71SMatthew G. Knepley 
250dce8aebaSBarry Smith   Fortran Note:
251dce8aebaSBarry Smith   From Fortran you must call `PetscQuadratureRestoreData()` when you are done with the data
2521fd49c25SBarry Smith 
253dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`, `PetscQuadratureSetData()`
25440d8ff71SMatthew G. Knepley @*/
255d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
256d71ae5a4SJacob Faibussowitsch {
25721454ff5SMatthew G. Knepley   PetscFunctionBegin;
2582cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
25921454ff5SMatthew G. Knepley   if (dim) {
260dadcf809SJacob Faibussowitsch     PetscValidIntPointer(dim, 2);
26121454ff5SMatthew G. Knepley     *dim = q->dim;
26221454ff5SMatthew G. Knepley   }
263a6b92713SMatthew G. Knepley   if (Nc) {
264dadcf809SJacob Faibussowitsch     PetscValidIntPointer(Nc, 3);
265a6b92713SMatthew G. Knepley     *Nc = q->Nc;
266a6b92713SMatthew G. Knepley   }
26721454ff5SMatthew G. Knepley   if (npoints) {
268dadcf809SJacob Faibussowitsch     PetscValidIntPointer(npoints, 4);
26921454ff5SMatthew G. Knepley     *npoints = q->numPoints;
27021454ff5SMatthew G. Knepley   }
27121454ff5SMatthew G. Knepley   if (points) {
272a6b92713SMatthew G. Knepley     PetscValidPointer(points, 5);
27321454ff5SMatthew G. Knepley     *points = q->points;
27421454ff5SMatthew G. Knepley   }
27521454ff5SMatthew G. Knepley   if (weights) {
276a6b92713SMatthew G. Knepley     PetscValidPointer(weights, 6);
27721454ff5SMatthew G. Knepley     *weights = q->weights;
27821454ff5SMatthew G. Knepley   }
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28021454ff5SMatthew G. Knepley }
28121454ff5SMatthew G. Knepley 
2824f9ab2b4SJed Brown /*@
2834f9ab2b4SJed Brown   PetscQuadratureEqual - determine whether two quadratures are equivalent
2844f9ab2b4SJed Brown 
2854f9ab2b4SJed Brown   Input Parameters:
286dce8aebaSBarry Smith + A - A `PetscQuadrature` object
287dce8aebaSBarry Smith - B - Another `PetscQuadrature` object
2884f9ab2b4SJed Brown 
2894f9ab2b4SJed Brown   Output Parameters:
290dce8aebaSBarry Smith . equal - `PETSC_TRUE` if the quadratures are the same
2914f9ab2b4SJed Brown 
2924f9ab2b4SJed Brown   Level: intermediate
2934f9ab2b4SJed Brown 
294dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`
2954f9ab2b4SJed Brown @*/
296d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureEqual(PetscQuadrature A, PetscQuadrature B, PetscBool *equal)
297d71ae5a4SJacob Faibussowitsch {
2984f9ab2b4SJed Brown   PetscFunctionBegin;
2994f9ab2b4SJed Brown   PetscValidHeaderSpecific(A, PETSCQUADRATURE_CLASSID, 1);
3004f9ab2b4SJed Brown   PetscValidHeaderSpecific(B, PETSCQUADRATURE_CLASSID, 2);
3014f9ab2b4SJed Brown   PetscValidBoolPointer(equal, 3);
3024f9ab2b4SJed Brown   *equal = PETSC_FALSE;
3033ba16761SJacob Faibussowitsch   if (A->dim != B->dim || A->Nc != B->Nc || A->order != B->order || A->numPoints != B->numPoints) PetscFunctionReturn(PETSC_SUCCESS);
3044f9ab2b4SJed Brown   for (PetscInt i = 0; i < A->numPoints * A->dim; i++) {
3053ba16761SJacob Faibussowitsch     if (PetscAbsReal(A->points[i] - B->points[i]) > PETSC_SMALL) PetscFunctionReturn(PETSC_SUCCESS);
3064f9ab2b4SJed Brown   }
3074f9ab2b4SJed Brown   if (!A->weights && !B->weights) {
3084f9ab2b4SJed Brown     *equal = PETSC_TRUE;
3093ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
3104f9ab2b4SJed Brown   }
3114f9ab2b4SJed Brown   if (A->weights && B->weights) {
3124f9ab2b4SJed Brown     for (PetscInt i = 0; i < A->numPoints; i++) {
3133ba16761SJacob Faibussowitsch       if (PetscAbsReal(A->weights[i] - B->weights[i]) > PETSC_SMALL) PetscFunctionReturn(PETSC_SUCCESS);
3144f9ab2b4SJed Brown     }
3154f9ab2b4SJed Brown     *equal = PETSC_TRUE;
3164f9ab2b4SJed Brown   }
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3184f9ab2b4SJed Brown }
3194f9ab2b4SJed Brown 
320d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[])
321d71ae5a4SJacob Faibussowitsch {
322907761f8SToby Isaac   PetscScalar *Js, *Jinvs;
323907761f8SToby Isaac   PetscInt     i, j, k;
324907761f8SToby Isaac   PetscBLASInt bm, bn, info;
325907761f8SToby Isaac 
326907761f8SToby Isaac   PetscFunctionBegin;
3273ba16761SJacob Faibussowitsch   if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS);
3289566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(m, &bm));
3299566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &bn));
330907761f8SToby Isaac #if defined(PETSC_USE_COMPLEX)
3319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m * n, &Js, m * n, &Jinvs));
33228222859SToby Isaac   for (i = 0; i < m * n; i++) Js[i] = J[i];
333907761f8SToby Isaac #else
334907761f8SToby Isaac   Js    = (PetscReal *)J;
335907761f8SToby Isaac   Jinvs = Jinv;
336907761f8SToby Isaac #endif
337907761f8SToby Isaac   if (m == n) {
338907761f8SToby Isaac     PetscBLASInt *pivots;
339907761f8SToby Isaac     PetscScalar  *W;
340907761f8SToby Isaac 
3419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(m, &pivots, m, &W));
342907761f8SToby Isaac 
3439566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(Jinvs, Js, m * m));
344792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info));
34563a3b9bcSJacob Faibussowitsch     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscInt_FMT, (PetscInt)info);
346792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info));
34763a3b9bcSJacob Faibussowitsch     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscInt_FMT, (PetscInt)info);
3489566063dSJacob Faibussowitsch     PetscCall(PetscFree2(pivots, W));
349907761f8SToby Isaac   } else if (m < n) {
350907761f8SToby Isaac     PetscScalar  *JJT;
351907761f8SToby Isaac     PetscBLASInt *pivots;
352907761f8SToby Isaac     PetscScalar  *W;
353907761f8SToby Isaac 
3549566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m * m, &JJT));
3559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(m, &pivots, m, &W));
356907761f8SToby Isaac     for (i = 0; i < m; i++) {
357907761f8SToby Isaac       for (j = 0; j < m; j++) {
358907761f8SToby Isaac         PetscScalar val = 0.;
359907761f8SToby Isaac 
360907761f8SToby Isaac         for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k];
361907761f8SToby Isaac         JJT[i * m + j] = val;
362907761f8SToby Isaac       }
363907761f8SToby Isaac     }
364907761f8SToby Isaac 
365792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info));
36663a3b9bcSJacob Faibussowitsch     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscInt_FMT, (PetscInt)info);
367792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info));
36863a3b9bcSJacob Faibussowitsch     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscInt_FMT, (PetscInt)info);
369907761f8SToby Isaac     for (i = 0; i < n; i++) {
370907761f8SToby Isaac       for (j = 0; j < m; j++) {
371907761f8SToby Isaac         PetscScalar val = 0.;
372907761f8SToby Isaac 
373907761f8SToby Isaac         for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j];
374907761f8SToby Isaac         Jinvs[i * m + j] = val;
375907761f8SToby Isaac       }
376907761f8SToby Isaac     }
3779566063dSJacob Faibussowitsch     PetscCall(PetscFree2(pivots, W));
3789566063dSJacob Faibussowitsch     PetscCall(PetscFree(JJT));
379907761f8SToby Isaac   } else {
380907761f8SToby Isaac     PetscScalar  *JTJ;
381907761f8SToby Isaac     PetscBLASInt *pivots;
382907761f8SToby Isaac     PetscScalar  *W;
383907761f8SToby Isaac 
3849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n * n, &JTJ));
3859566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n, &pivots, n, &W));
386907761f8SToby Isaac     for (i = 0; i < n; i++) {
387907761f8SToby Isaac       for (j = 0; j < n; j++) {
388907761f8SToby Isaac         PetscScalar val = 0.;
389907761f8SToby Isaac 
390907761f8SToby Isaac         for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j];
391907761f8SToby Isaac         JTJ[i * n + j] = val;
392907761f8SToby Isaac       }
393907761f8SToby Isaac     }
394907761f8SToby Isaac 
395792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bn, pivots, &info));
39663a3b9bcSJacob Faibussowitsch     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscInt_FMT, (PetscInt)info);
397792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info));
39863a3b9bcSJacob Faibussowitsch     PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscInt_FMT, (PetscInt)info);
399907761f8SToby Isaac     for (i = 0; i < n; i++) {
400907761f8SToby Isaac       for (j = 0; j < m; j++) {
401907761f8SToby Isaac         PetscScalar val = 0.;
402907761f8SToby Isaac 
403907761f8SToby Isaac         for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k];
404907761f8SToby Isaac         Jinvs[i * m + j] = val;
405907761f8SToby Isaac       }
406907761f8SToby Isaac     }
4079566063dSJacob Faibussowitsch     PetscCall(PetscFree2(pivots, W));
4089566063dSJacob Faibussowitsch     PetscCall(PetscFree(JTJ));
409907761f8SToby Isaac   }
410907761f8SToby Isaac #if defined(PETSC_USE_COMPLEX)
41128222859SToby Isaac   for (i = 0; i < m * n; i++) Jinv[i] = PetscRealPart(Jinvs[i]);
4129566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Js, Jinvs));
413907761f8SToby Isaac #endif
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
415907761f8SToby Isaac }
416907761f8SToby Isaac 
417907761f8SToby Isaac /*@
418907761f8SToby Isaac    PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation.
419907761f8SToby Isaac 
42020f4b53cSBarry Smith    Collective
421907761f8SToby Isaac 
4224165533cSJose E. Roman    Input Parameters:
423907761f8SToby Isaac +  q - the quadrature functional
424907761f8SToby Isaac .  imageDim - the dimension of the image of the transformation
425907761f8SToby Isaac .  origin - a point in the original space
426907761f8SToby Isaac .  originImage - the image of the origin under the transformation
427907761f8SToby Isaac .  J - the Jacobian of the image: an [imageDim x dim] matrix in row major order
428dce8aebaSBarry Smith -  formDegree - transform the quadrature weights as k-forms of this form degree (if the number of components is a multiple of (dim choose formDegree), it is assumed that they represent multiple k-forms) [see `PetscDTAltVPullback()` for interpretation of formDegree]
429907761f8SToby Isaac 
4304165533cSJose E. Roman    Output Parameters:
431907761f8SToby Isaac .  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 been pulled-back by the pseudoinverse of J to the k-form weights in the image space.
432907761f8SToby Isaac 
4336c877ef6SSatish Balay    Level: intermediate
4346c877ef6SSatish Balay 
435dce8aebaSBarry Smith    Note:
436dce8aebaSBarry Smith    The new quadrature rule will have a different number of components if spaces have different dimensions.  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.
437dce8aebaSBarry Smith 
438dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
439907761f8SToby Isaac @*/
440d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formDegree, PetscQuadrature *Jinvstarq)
441d71ae5a4SJacob Faibussowitsch {
442907761f8SToby Isaac   PetscInt         dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c;
443907761f8SToby Isaac   const PetscReal *points;
444907761f8SToby Isaac   const PetscReal *weights;
445907761f8SToby Isaac   PetscReal       *imagePoints, *imageWeights;
446907761f8SToby Isaac   PetscReal       *Jinv;
447907761f8SToby Isaac   PetscReal       *Jinvstar;
448907761f8SToby Isaac 
449907761f8SToby Isaac   PetscFunctionBegin;
450d4afb720SToby Isaac   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
45163a3b9bcSJacob Faibussowitsch   PetscCheck(imageDim >= PetscAbsInt(formDegree), PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %" PetscInt_FMT "-form in %" PetscInt_FMT " dimensions", PetscAbsInt(formDegree), imageDim);
4529566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights));
4539566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &formSize));
45463a3b9bcSJacob Faibussowitsch   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);
455907761f8SToby Isaac   Ncopies = Nc / formSize;
4569566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(imageDim, PetscAbsInt(formDegree), &imageFormSize));
457907761f8SToby Isaac   imageNc = Ncopies * imageFormSize;
4589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Npoints * imageDim, &imagePoints));
4599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Npoints * imageNc, &imageWeights));
4609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar));
4619566063dSJacob Faibussowitsch   PetscCall(PetscDTJacobianInverse_Internal(imageDim, dim, J, Jinv));
4629566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formDegree, Jinvstar));
463907761f8SToby Isaac   for (pt = 0; pt < Npoints; pt++) {
464907761f8SToby Isaac     const PetscReal *point      = &points[pt * dim];
465907761f8SToby Isaac     PetscReal       *imagePoint = &imagePoints[pt * imageDim];
466907761f8SToby Isaac 
467907761f8SToby Isaac     for (i = 0; i < imageDim; i++) {
468907761f8SToby Isaac       PetscReal val = originImage[i];
469907761f8SToby Isaac 
470907761f8SToby Isaac       for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]);
471907761f8SToby Isaac       imagePoint[i] = val;
472907761f8SToby Isaac     }
473907761f8SToby Isaac     for (c = 0; c < Ncopies; c++) {
474907761f8SToby Isaac       const PetscReal *form      = &weights[pt * Nc + c * formSize];
475907761f8SToby Isaac       PetscReal       *imageForm = &imageWeights[pt * imageNc + c * imageFormSize];
476907761f8SToby Isaac 
477907761f8SToby Isaac       for (i = 0; i < imageFormSize; i++) {
478907761f8SToby Isaac         PetscReal val = 0.;
479907761f8SToby Isaac 
480907761f8SToby Isaac         for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j];
481907761f8SToby Isaac         imageForm[i] = val;
482907761f8SToby Isaac       }
483907761f8SToby Isaac     }
484907761f8SToby Isaac   }
4859566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq));
4869566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights));
4879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Jinv, Jinvstar));
4883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
489907761f8SToby Isaac }
490907761f8SToby Isaac 
49140d8ff71SMatthew G. Knepley /*@C
49240d8ff71SMatthew G. Knepley   PetscQuadratureSetData - Sets the data defining the quadrature
49340d8ff71SMatthew G. Knepley 
49420f4b53cSBarry Smith   Not Collective
49540d8ff71SMatthew G. Knepley 
49640d8ff71SMatthew G. Knepley   Input Parameters:
497dce8aebaSBarry Smith + q  - The `PetscQuadrature` object
49840d8ff71SMatthew G. Knepley . dim - The spatial dimension
499e2b35d93SBarry Smith . Nc - The number of components
50040d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
50140d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
50240d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
50340d8ff71SMatthew G. Knepley 
50440d8ff71SMatthew G. Knepley   Level: intermediate
50540d8ff71SMatthew G. Knepley 
506dce8aebaSBarry Smith   Note:
507dce8aebaSBarry Smith   This routine owns the references to points and weights, so they must be allocated using `PetscMalloc()` and the user should not free them.
508dce8aebaSBarry Smith 
509dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscQuadratureCreate()`, `PetscQuadratureGetData()`
51040d8ff71SMatthew G. Knepley @*/
511d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
512d71ae5a4SJacob Faibussowitsch {
51321454ff5SMatthew G. Knepley   PetscFunctionBegin;
5142cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
51521454ff5SMatthew G. Knepley   if (dim >= 0) q->dim = dim;
516a6b92713SMatthew G. Knepley   if (Nc >= 0) q->Nc = Nc;
51721454ff5SMatthew G. Knepley   if (npoints >= 0) q->numPoints = npoints;
51821454ff5SMatthew G. Knepley   if (points) {
519dadcf809SJacob Faibussowitsch     PetscValidRealPointer(points, 5);
52021454ff5SMatthew G. Knepley     q->points = points;
52121454ff5SMatthew G. Knepley   }
52221454ff5SMatthew G. Knepley   if (weights) {
523dadcf809SJacob Faibussowitsch     PetscValidRealPointer(weights, 6);
52421454ff5SMatthew G. Knepley     q->weights = weights;
52521454ff5SMatthew G. Knepley   }
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
527f9fd7fdbSMatthew G. Knepley }
528f9fd7fdbSMatthew G. Knepley 
529d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v)
530d71ae5a4SJacob Faibussowitsch {
531d9bac1caSLisandro Dalcin   PetscInt          q, d, c;
532d9bac1caSLisandro Dalcin   PetscViewerFormat format;
533d9bac1caSLisandro Dalcin 
534d9bac1caSLisandro Dalcin   PetscFunctionBegin;
53563a3b9bcSJacob Faibussowitsch   if (quad->Nc > 1) PetscCall(PetscViewerASCIIPrintf(v, "Quadrature of order %" PetscInt_FMT " on %" PetscInt_FMT " points (dim %" PetscInt_FMT ") with %" PetscInt_FMT " components\n", quad->order, quad->numPoints, quad->dim, quad->Nc));
53663a3b9bcSJacob Faibussowitsch   else PetscCall(PetscViewerASCIIPrintf(v, "Quadrature of order %" PetscInt_FMT " on %" PetscInt_FMT " points (dim %" PetscInt_FMT ")\n", quad->order, quad->numPoints, quad->dim));
5379566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(v, &format));
5383ba16761SJacob Faibussowitsch   if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS);
539d9bac1caSLisandro Dalcin   for (q = 0; q < quad->numPoints; ++q) {
54063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "p%" PetscInt_FMT " (", q));
5419566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(v, PETSC_FALSE));
542d9bac1caSLisandro Dalcin     for (d = 0; d < quad->dim; ++d) {
5439566063dSJacob Faibussowitsch       if (d) PetscCall(PetscViewerASCIIPrintf(v, ", "));
5449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q * quad->dim + d]));
545d9bac1caSLisandro Dalcin     }
5469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, ") "));
54763a3b9bcSJacob Faibussowitsch     if (quad->Nc > 1) PetscCall(PetscViewerASCIIPrintf(v, "w%" PetscInt_FMT " (", q));
548d9bac1caSLisandro Dalcin     for (c = 0; c < quad->Nc; ++c) {
5499566063dSJacob Faibussowitsch       if (c) PetscCall(PetscViewerASCIIPrintf(v, ", "));
5509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q * quad->Nc + c]));
551d9bac1caSLisandro Dalcin     }
5529566063dSJacob Faibussowitsch     if (quad->Nc > 1) PetscCall(PetscViewerASCIIPrintf(v, ")"));
5539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "\n"));
5549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(v, PETSC_TRUE));
555d9bac1caSLisandro Dalcin   }
5563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
557d9bac1caSLisandro Dalcin }
558d9bac1caSLisandro Dalcin 
55940d8ff71SMatthew G. Knepley /*@C
560dce8aebaSBarry Smith   PetscQuadratureView - View a `PetscQuadrature` object
56140d8ff71SMatthew G. Knepley 
56220f4b53cSBarry Smith   Collective
56340d8ff71SMatthew G. Knepley 
56440d8ff71SMatthew G. Knepley   Input Parameters:
565dce8aebaSBarry Smith + quad  - The `PetscQuadrature` object
566dce8aebaSBarry Smith - viewer - The `PetscViewer` object
56740d8ff71SMatthew G. Knepley 
56840d8ff71SMatthew G. Knepley   Level: beginner
56940d8ff71SMatthew G. Knepley 
570dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscViewer`, `PetscQuadratureCreate()`, `PetscQuadratureGetData()`
57140d8ff71SMatthew G. Knepley @*/
572d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
573d71ae5a4SJacob Faibussowitsch {
574d9bac1caSLisandro Dalcin   PetscBool iascii;
575f9fd7fdbSMatthew G. Knepley 
576f9fd7fdbSMatthew G. Knepley   PetscFunctionBegin;
577d9bac1caSLisandro Dalcin   PetscValidHeader(quad, 1);
578d9bac1caSLisandro Dalcin   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
5799566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)quad), &viewer));
5809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
5819566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
5829566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscQuadratureView_Ascii(quad, viewer));
5839566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
5843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
585bfa639d9SMatthew G. Knepley }
586bfa639d9SMatthew G. Knepley 
58789710940SMatthew G. Knepley /*@C
58889710940SMatthew G. Knepley   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
58989710940SMatthew G. Knepley 
59020f4b53cSBarry Smith   Not Collective; No Fortran Support
59189710940SMatthew G. Knepley 
592d8d19677SJose E. Roman   Input Parameters:
593dce8aebaSBarry Smith + q - The original `PetscQuadrature`
59489710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into
59589710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement
59689710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement
59789710940SMatthew G. Knepley 
59889710940SMatthew G. Knepley   Output Parameters:
59989710940SMatthew G. Knepley . dim - The dimension
60089710940SMatthew G. Knepley 
60120f4b53cSBarry Smith   Level: intermediate
60220f4b53cSBarry Smith 
603dce8aebaSBarry Smith   Note:
604dce8aebaSBarry Smith   Together v0 and jac define an affine mapping from the original reference element to each subelement
60589710940SMatthew G. Knepley 
606dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscFECreate()`, `PetscSpaceGetDimension()`, `PetscDualSpaceGetDimension()`
60789710940SMatthew G. Knepley @*/
608d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
609d71ae5a4SJacob Faibussowitsch {
61089710940SMatthew G. Knepley   const PetscReal *points, *weights;
61189710940SMatthew G. Knepley   PetscReal       *pointsRef, *weightsRef;
612a6b92713SMatthew G. Knepley   PetscInt         dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
61389710940SMatthew G. Knepley 
61489710940SMatthew G. Knepley   PetscFunctionBegin;
6152cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
616dadcf809SJacob Faibussowitsch   PetscValidRealPointer(v0, 3);
617dadcf809SJacob Faibussowitsch   PetscValidRealPointer(jac, 4);
61889710940SMatthew G. Knepley   PetscValidPointer(qref, 5);
6199566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, qref));
6209566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetOrder(q, &order));
6219566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights));
62289710940SMatthew G. Knepley   npointsRef = npoints * numSubelements;
6239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npointsRef * dim, &pointsRef));
6249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npointsRef * Nc, &weightsRef));
62589710940SMatthew G. Knepley   for (c = 0; c < numSubelements; ++c) {
62689710940SMatthew G. Knepley     for (p = 0; p < npoints; ++p) {
62789710940SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
62889710940SMatthew G. Knepley         pointsRef[(c * npoints + p) * dim + d] = v0[c * dim + d];
629ad540459SPierre Jolivet         for (e = 0; e < dim; ++e) pointsRef[(c * npoints + p) * dim + d] += jac[(c * dim + d) * dim + e] * (points[p * dim + e] + 1.0);
63089710940SMatthew G. Knepley       }
63189710940SMatthew G. Knepley       /* Could also use detJ here */
632a6b92713SMatthew G. Knepley       for (cp = 0; cp < Nc; ++cp) weightsRef[(c * npoints + p) * Nc + cp] = weights[p * Nc + cp] / numSubelements;
63389710940SMatthew G. Knepley     }
63489710940SMatthew G. Knepley   }
6359566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetOrder(*qref, order));
6369566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef));
6373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63889710940SMatthew G. Knepley }
63989710940SMatthew G. Knepley 
64094e21283SToby Isaac /* Compute the coefficients for the Jacobi polynomial recurrence,
64194e21283SToby Isaac  *
64294e21283SToby Isaac  * J^{a,b}_n(x) = (cnm1 + cnm1x * x) * J^{a,b}_{n-1}(x) - cnm2 * J^{a,b}_{n-2}(x).
64394e21283SToby Isaac  */
64494e21283SToby Isaac #define PetscDTJacobiRecurrence_Internal(n, a, b, cnm1, cnm1x, cnm2) \
64594e21283SToby Isaac   do { \
64694e21283SToby Isaac     PetscReal _a = (a); \
64794e21283SToby Isaac     PetscReal _b = (b); \
64894e21283SToby Isaac     PetscReal _n = (n); \
64994e21283SToby Isaac     if (n == 1) { \
65094e21283SToby Isaac       (cnm1)  = (_a - _b) * 0.5; \
65194e21283SToby Isaac       (cnm1x) = (_a + _b + 2.) * 0.5; \
65294e21283SToby Isaac       (cnm2)  = 0.; \
65394e21283SToby Isaac     } else { \
65494e21283SToby Isaac       PetscReal _2n  = _n + _n; \
65594e21283SToby Isaac       PetscReal _d   = (_2n * (_n + _a + _b) * (_2n + _a + _b - 2)); \
65694e21283SToby Isaac       PetscReal _n1  = (_2n + _a + _b - 1.) * (_a * _a - _b * _b); \
65794e21283SToby Isaac       PetscReal _n1x = (_2n + _a + _b - 1.) * (_2n + _a + _b) * (_2n + _a + _b - 2); \
65894e21283SToby Isaac       PetscReal _n2  = 2. * ((_n + _a - 1.) * (_n + _b - 1.) * (_2n + _a + _b)); \
65994e21283SToby Isaac       (cnm1)         = _n1 / _d; \
66094e21283SToby Isaac       (cnm1x)        = _n1x / _d; \
66194e21283SToby Isaac       (cnm2)         = _n2 / _d; \
66294e21283SToby Isaac     } \
66394e21283SToby Isaac   } while (0)
66494e21283SToby Isaac 
665fbdc3dfeSToby Isaac /*@
666fbdc3dfeSToby Isaac   PetscDTJacobiNorm - Compute the weighted L2 norm of a Jacobi polynomial.
667fbdc3dfeSToby Isaac 
668fbdc3dfeSToby Isaac   $\| P^{\alpha,\beta}_n \|_{\alpha,\beta}^2 = \int_{-1}^1 (1 + x)^{\alpha} (1 - x)^{\beta} P^{\alpha,\beta}_n (x)^2 dx.$
669fbdc3dfeSToby Isaac 
6704165533cSJose E. Roman   Input Parameters:
671fbdc3dfeSToby Isaac - alpha - the left exponent > -1
672fbdc3dfeSToby Isaac . beta - the right exponent > -1
673fbdc3dfeSToby Isaac + n - the polynomial degree
674fbdc3dfeSToby Isaac 
6754165533cSJose E. Roman   Output Parameter:
676fbdc3dfeSToby Isaac . norm - the weighted L2 norm
677fbdc3dfeSToby Isaac 
678fbdc3dfeSToby Isaac   Level: beginner
679fbdc3dfeSToby Isaac 
680dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscDTJacobiEval()`
681fbdc3dfeSToby Isaac @*/
682d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTJacobiNorm(PetscReal alpha, PetscReal beta, PetscInt n, PetscReal *norm)
683d71ae5a4SJacob Faibussowitsch {
684fbdc3dfeSToby Isaac   PetscReal twoab1;
685fbdc3dfeSToby Isaac   PetscReal gr;
686fbdc3dfeSToby Isaac 
687fbdc3dfeSToby Isaac   PetscFunctionBegin;
68808401ef6SPierre Jolivet   PetscCheck(alpha > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exponent alpha %g <= -1. invalid", (double)alpha);
68908401ef6SPierre Jolivet   PetscCheck(beta > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exponent beta %g <= -1. invalid", (double)beta);
69063a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "n %" PetscInt_FMT " < 0 invalid", n);
691fbdc3dfeSToby Isaac   twoab1 = PetscPowReal(2., alpha + beta + 1.);
692fbdc3dfeSToby Isaac #if defined(PETSC_HAVE_LGAMMA)
693fbdc3dfeSToby Isaac   if (!n) {
694fbdc3dfeSToby Isaac     gr = PetscExpReal(PetscLGamma(alpha + 1.) + PetscLGamma(beta + 1.) - PetscLGamma(alpha + beta + 2.));
695fbdc3dfeSToby Isaac   } else {
696fbdc3dfeSToby Isaac     gr = PetscExpReal(PetscLGamma(n + alpha + 1.) + PetscLGamma(n + beta + 1.) - (PetscLGamma(n + 1.) + PetscLGamma(n + alpha + beta + 1.))) / (n + n + alpha + beta + 1.);
697fbdc3dfeSToby Isaac   }
698fbdc3dfeSToby Isaac #else
699fbdc3dfeSToby Isaac   {
700fbdc3dfeSToby Isaac     PetscInt alphai = (PetscInt)alpha;
701fbdc3dfeSToby Isaac     PetscInt betai  = (PetscInt)beta;
702fbdc3dfeSToby Isaac     PetscInt i;
703fbdc3dfeSToby Isaac 
704fbdc3dfeSToby Isaac     gr = n ? (1. / (n + n + alpha + beta + 1.)) : 1.;
705fbdc3dfeSToby Isaac     if ((PetscReal)alphai == alpha) {
706fbdc3dfeSToby Isaac       if (!n) {
707fbdc3dfeSToby Isaac         for (i = 0; i < alphai; i++) gr *= (i + 1.) / (beta + i + 1.);
708fbdc3dfeSToby Isaac         gr /= (alpha + beta + 1.);
709fbdc3dfeSToby Isaac       } else {
710fbdc3dfeSToby Isaac         for (i = 0; i < alphai; i++) gr *= (n + i + 1.) / (n + beta + i + 1.);
711fbdc3dfeSToby Isaac       }
712fbdc3dfeSToby Isaac     } else if ((PetscReal)betai == beta) {
713fbdc3dfeSToby Isaac       if (!n) {
714fbdc3dfeSToby Isaac         for (i = 0; i < betai; i++) gr *= (i + 1.) / (alpha + i + 2.);
715fbdc3dfeSToby Isaac         gr /= (alpha + beta + 1.);
716fbdc3dfeSToby Isaac       } else {
717fbdc3dfeSToby Isaac         for (i = 0; i < betai; i++) gr *= (n + i + 1.) / (n + alpha + i + 1.);
718fbdc3dfeSToby Isaac       }
719fbdc3dfeSToby Isaac     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "lgamma() - math routine is unavailable.");
720fbdc3dfeSToby Isaac   }
721fbdc3dfeSToby Isaac #endif
722fbdc3dfeSToby Isaac   *norm = PetscSqrtReal(twoab1 * gr);
7233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
724fbdc3dfeSToby Isaac }
725fbdc3dfeSToby Isaac 
726d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTJacobiEval_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscInt k, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *p)
727d71ae5a4SJacob Faibussowitsch {
72894e21283SToby Isaac   PetscReal ak, bk;
72994e21283SToby Isaac   PetscReal abk1;
73094e21283SToby Isaac   PetscInt  i, l, maxdegree;
73194e21283SToby Isaac 
73294e21283SToby Isaac   PetscFunctionBegin;
73394e21283SToby Isaac   maxdegree = degrees[ndegree - 1] - k;
73494e21283SToby Isaac   ak        = a + k;
73594e21283SToby Isaac   bk        = b + k;
73694e21283SToby Isaac   abk1      = a + b + k + 1.;
73794e21283SToby Isaac   if (maxdegree < 0) {
7389371c9d4SSatish Balay     for (i = 0; i < npoints; i++)
7399371c9d4SSatish Balay       for (l = 0; l < ndegree; l++) p[i * ndegree + l] = 0.;
7403ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
74194e21283SToby Isaac   }
74294e21283SToby Isaac   for (i = 0; i < npoints; i++) {
74394e21283SToby Isaac     PetscReal pm1, pm2, x;
74494e21283SToby Isaac     PetscReal cnm1, cnm1x, cnm2;
74594e21283SToby Isaac     PetscInt  j, m;
74694e21283SToby Isaac 
74794e21283SToby Isaac     x   = points[i];
74894e21283SToby Isaac     pm2 = 1.;
74994e21283SToby Isaac     PetscDTJacobiRecurrence_Internal(1, ak, bk, cnm1, cnm1x, cnm2);
75094e21283SToby Isaac     pm1 = (cnm1 + cnm1x * x);
75194e21283SToby Isaac     l   = 0;
752ad540459SPierre Jolivet     while (l < ndegree && degrees[l] - k < 0) p[l++] = 0.;
75394e21283SToby Isaac     while (l < ndegree && degrees[l] - k == 0) {
75494e21283SToby Isaac       p[l] = pm2;
75594e21283SToby Isaac       for (m = 0; m < k; m++) p[l] *= (abk1 + m) * 0.5;
75694e21283SToby Isaac       l++;
75794e21283SToby Isaac     }
75894e21283SToby Isaac     while (l < ndegree && degrees[l] - k == 1) {
75994e21283SToby Isaac       p[l] = pm1;
76094e21283SToby Isaac       for (m = 0; m < k; m++) p[l] *= (abk1 + 1 + m) * 0.5;
76194e21283SToby Isaac       l++;
76294e21283SToby Isaac     }
76394e21283SToby Isaac     for (j = 2; j <= maxdegree; j++) {
76494e21283SToby Isaac       PetscReal pp;
76594e21283SToby Isaac 
76694e21283SToby Isaac       PetscDTJacobiRecurrence_Internal(j, ak, bk, cnm1, cnm1x, cnm2);
76794e21283SToby Isaac       pp  = (cnm1 + cnm1x * x) * pm1 - cnm2 * pm2;
76894e21283SToby Isaac       pm2 = pm1;
76994e21283SToby Isaac       pm1 = pp;
77094e21283SToby Isaac       while (l < ndegree && degrees[l] - k == j) {
77194e21283SToby Isaac         p[l] = pp;
77294e21283SToby Isaac         for (m = 0; m < k; m++) p[l] *= (abk1 + j + m) * 0.5;
77394e21283SToby Isaac         l++;
77494e21283SToby Isaac       }
77594e21283SToby Isaac     }
77694e21283SToby Isaac     p += ndegree;
77794e21283SToby Isaac   }
7783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
77994e21283SToby Isaac }
78094e21283SToby Isaac 
78137045ce4SJed Brown /*@
782dce8aebaSBarry Smith   PetscDTJacobiEvalJet - Evaluate the jet (function and derivatives) of the Jacobi polynomials basis up to a given degree.
783dce8aebaSBarry Smith   The Jacobi polynomials with indices $\alpha$ and $\beta$ are orthogonal with respect to the weighted inner product
784dce8aebaSBarry Smith   $\langle f, g \rangle = \int_{-1}^1 (1+x)^{\alpha} (1-x)^{\beta} f(x) g(x) dx$.
785fbdc3dfeSToby Isaac 
7864165533cSJose E. Roman   Input Parameters:
787fbdc3dfeSToby Isaac + alpha - the left exponent of the weight
788fbdc3dfeSToby Isaac . beta - the right exponetn of the weight
789fbdc3dfeSToby Isaac . npoints - the number of points to evaluate the polynomials at
790fbdc3dfeSToby Isaac . points - [npoints] array of point coordinates
791fbdc3dfeSToby Isaac . degree - the maximm degree polynomial space to evaluate, (degree + 1) will be evaluated total.
792fbdc3dfeSToby Isaac - k - the maximum derivative to evaluate in the jet, (k + 1) will be evaluated total.
793fbdc3dfeSToby Isaac 
7946aad120cSJose E. Roman   Output Parameters:
795fbdc3dfeSToby Isaac - p - an array containing the evaluations of the Jacobi polynomials's jets on the points.  the size is (degree + 1) x
796fbdc3dfeSToby Isaac   (k + 1) x npoints, which also describes the order of the dimensions of this three-dimensional array: the first
797fbdc3dfeSToby Isaac   (slowest varying) dimension is polynomial degree; the second dimension is derivative order; the third (fastest
798fbdc3dfeSToby Isaac   varying) dimension is the index of the evaluation point.
799fbdc3dfeSToby Isaac 
800fbdc3dfeSToby Isaac   Level: advanced
801fbdc3dfeSToby Isaac 
802db781477SPatrick Sanan .seealso: `PetscDTJacobiEval()`, `PetscDTPKDEvalJet()`
803fbdc3dfeSToby Isaac @*/
804d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTJacobiEvalJet(PetscReal alpha, PetscReal beta, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt k, PetscReal p[])
805d71ae5a4SJacob Faibussowitsch {
806fbdc3dfeSToby Isaac   PetscInt   i, j, l;
807fbdc3dfeSToby Isaac   PetscInt  *degrees;
808fbdc3dfeSToby Isaac   PetscReal *psingle;
809fbdc3dfeSToby Isaac 
810fbdc3dfeSToby Isaac   PetscFunctionBegin;
811fbdc3dfeSToby Isaac   if (degree == 0) {
812fbdc3dfeSToby Isaac     PetscInt zero = 0;
813fbdc3dfeSToby Isaac 
81448a46eb9SPierre Jolivet     for (i = 0; i <= k; i++) PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, i, points, 1, &zero, &p[i * npoints]));
8153ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
816fbdc3dfeSToby Isaac   }
8179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(degree + 1, &degrees));
8189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((degree + 1) * npoints, &psingle));
819fbdc3dfeSToby Isaac   for (i = 0; i <= degree; i++) degrees[i] = i;
820fbdc3dfeSToby Isaac   for (i = 0; i <= k; i++) {
8219566063dSJacob Faibussowitsch     PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, i, points, degree + 1, degrees, psingle));
822fbdc3dfeSToby Isaac     for (j = 0; j <= degree; j++) {
823ad540459SPierre Jolivet       for (l = 0; l < npoints; l++) p[(j * (k + 1) + i) * npoints + l] = psingle[l * (degree + 1) + j];
824fbdc3dfeSToby Isaac     }
825fbdc3dfeSToby Isaac   }
8269566063dSJacob Faibussowitsch   PetscCall(PetscFree(psingle));
8279566063dSJacob Faibussowitsch   PetscCall(PetscFree(degrees));
8283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
829fbdc3dfeSToby Isaac }
830fbdc3dfeSToby Isaac 
831fbdc3dfeSToby Isaac /*@
832dce8aebaSBarry Smith    PetscDTJacobiEval - evaluate Jacobi polynomials for the weight function $(1.+x)^{\alpha} (1.-x)^{\beta}$ at a set of points
83394e21283SToby Isaac                        at points
83494e21283SToby Isaac 
83594e21283SToby Isaac    Not Collective
83694e21283SToby Isaac 
8374165533cSJose E. Roman    Input Parameters:
83894e21283SToby Isaac +  npoints - number of spatial points to evaluate at
83994e21283SToby Isaac .  alpha - the left exponent > -1
84094e21283SToby Isaac .  beta - the right exponent > -1
84194e21283SToby Isaac .  points - array of locations to evaluate at
84294e21283SToby Isaac .  ndegree - number of basis degrees to evaluate
84394e21283SToby Isaac -  degrees - sorted array of degrees to evaluate
84494e21283SToby Isaac 
8454165533cSJose E. Roman    Output Parameters:
84694e21283SToby Isaac +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
84794e21283SToby Isaac .  D - row-oriented derivative evaluation matrix (or NULL)
84894e21283SToby Isaac -  D2 - row-oriented second derivative evaluation matrix (or NULL)
84994e21283SToby Isaac 
85094e21283SToby Isaac    Level: intermediate
85194e21283SToby Isaac 
852dce8aebaSBarry Smith .seealso: `PetscDTGaussQuadrature()`, `PetscDTLegendreEval()`
85394e21283SToby Isaac @*/
854d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTJacobiEval(PetscInt npoints, PetscReal alpha, PetscReal beta, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *B, PetscReal *D, PetscReal *D2)
855d71ae5a4SJacob Faibussowitsch {
85694e21283SToby Isaac   PetscFunctionBegin;
85708401ef6SPierre Jolivet   PetscCheck(alpha > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "alpha must be > -1.");
85808401ef6SPierre Jolivet   PetscCheck(beta > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "beta must be > -1.");
8593ba16761SJacob Faibussowitsch   if (!npoints || !ndegree) PetscFunctionReturn(PETSC_SUCCESS);
8609566063dSJacob Faibussowitsch   if (B) PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, 0, points, ndegree, degrees, B));
8619566063dSJacob Faibussowitsch   if (D) PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, 1, points, ndegree, degrees, D));
8629566063dSJacob Faibussowitsch   if (D2) PetscCall(PetscDTJacobiEval_Internal(npoints, alpha, beta, 2, points, ndegree, degrees, D2));
8633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86494e21283SToby Isaac }
86594e21283SToby Isaac 
86694e21283SToby Isaac /*@
86794e21283SToby Isaac    PetscDTLegendreEval - evaluate Legendre polynomials at points
86837045ce4SJed Brown 
86937045ce4SJed Brown    Not Collective
87037045ce4SJed Brown 
8714165533cSJose E. Roman    Input Parameters:
87237045ce4SJed Brown +  npoints - number of spatial points to evaluate at
87337045ce4SJed Brown .  points - array of locations to evaluate at
87437045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
87537045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
87637045ce4SJed Brown 
8774165533cSJose E. Roman    Output Parameters:
8780298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
8790298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
8800298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
88137045ce4SJed Brown 
88237045ce4SJed Brown    Level: intermediate
88337045ce4SJed Brown 
884db781477SPatrick Sanan .seealso: `PetscDTGaussQuadrature()`
88537045ce4SJed Brown @*/
886d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTLegendreEval(PetscInt npoints, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *B, PetscReal *D, PetscReal *D2)
887d71ae5a4SJacob Faibussowitsch {
88837045ce4SJed Brown   PetscFunctionBegin;
8899566063dSJacob Faibussowitsch   PetscCall(PetscDTJacobiEval(npoints, 0., 0., points, ndegree, degrees, B, D, D2));
8903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
89137045ce4SJed Brown }
89237045ce4SJed Brown 
893fbdc3dfeSToby Isaac /*@
894fbdc3dfeSToby Isaac   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, then the index of x is smaller than the index of y)
895fbdc3dfeSToby Isaac 
896fbdc3dfeSToby Isaac   Input Parameters:
897fbdc3dfeSToby Isaac + len - the desired length of the degree tuple
898fbdc3dfeSToby Isaac - index - the index to convert: should be >= 0
899fbdc3dfeSToby Isaac 
900fbdc3dfeSToby Isaac   Output Parameter:
901fbdc3dfeSToby Isaac . degtup - will be filled with a tuple of degrees
902fbdc3dfeSToby Isaac 
903fbdc3dfeSToby Isaac   Level: beginner
904fbdc3dfeSToby Isaac 
905dce8aebaSBarry Smith   Note:
906dce8aebaSBarry Smith   For two tuples x and y with the same degree sum, partial degree sums over the final elements of the tuples
907fbdc3dfeSToby Isaac   acts as a tiebreaker.  For example, (2, 1, 1) and (1, 2, 1) have the same degree sum, but the degree sum over the
908fbdc3dfeSToby Isaac   last two elements is smaller for the former, so (2, 1, 1) < (1, 2, 1).
909fbdc3dfeSToby Isaac 
910db781477SPatrick Sanan .seealso: `PetscDTGradedOrderToIndex()`
911fbdc3dfeSToby Isaac @*/
912d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTIndexToGradedOrder(PetscInt len, PetscInt index, PetscInt degtup[])
913d71ae5a4SJacob Faibussowitsch {
914fbdc3dfeSToby Isaac   PetscInt i, total;
915fbdc3dfeSToby Isaac   PetscInt sum;
916fbdc3dfeSToby Isaac 
917fbdc3dfeSToby Isaac   PetscFunctionBeginHot;
91808401ef6SPierre Jolivet   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
91908401ef6SPierre Jolivet   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index must be non-negative");
920fbdc3dfeSToby Isaac   total = 1;
921fbdc3dfeSToby Isaac   sum   = 0;
922fbdc3dfeSToby Isaac   while (index >= total) {
923fbdc3dfeSToby Isaac     index -= total;
924fbdc3dfeSToby Isaac     total = (total * (len + sum)) / (sum + 1);
925fbdc3dfeSToby Isaac     sum++;
926fbdc3dfeSToby Isaac   }
927fbdc3dfeSToby Isaac   for (i = 0; i < len; i++) {
928fbdc3dfeSToby Isaac     PetscInt c;
929fbdc3dfeSToby Isaac 
930fbdc3dfeSToby Isaac     degtup[i] = sum;
931fbdc3dfeSToby Isaac     for (c = 0, total = 1; c < sum; c++) {
932fbdc3dfeSToby Isaac       /* going into the loop, total is the number of way to have a tuple of sum exactly c with length len - 1 - i */
933fbdc3dfeSToby Isaac       if (index < total) break;
934fbdc3dfeSToby Isaac       index -= total;
935fbdc3dfeSToby Isaac       total = (total * (len - 1 - i + c)) / (c + 1);
936fbdc3dfeSToby Isaac       degtup[i]--;
937fbdc3dfeSToby Isaac     }
938fbdc3dfeSToby Isaac     sum -= degtup[i];
939fbdc3dfeSToby Isaac   }
9403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
941fbdc3dfeSToby Isaac }
942fbdc3dfeSToby Isaac 
943fbdc3dfeSToby Isaac /*@
944dce8aebaSBarry Smith   PetscDTGradedOrderToIndex - convert a tuple into an index in a graded order, the inverse of `PetscDTIndexToGradedOrder()`.
945fbdc3dfeSToby Isaac 
946fbdc3dfeSToby Isaac   Input Parameters:
947fbdc3dfeSToby Isaac + len - the length of the degree tuple
948fbdc3dfeSToby Isaac - degtup - tuple with this length
949fbdc3dfeSToby Isaac 
950fbdc3dfeSToby Isaac   Output Parameter:
951fbdc3dfeSToby Isaac . index - index in graded order: >= 0
952fbdc3dfeSToby Isaac 
953fbdc3dfeSToby Isaac   Level: Beginner
954fbdc3dfeSToby Isaac 
955dce8aebaSBarry Smith   Note:
956dce8aebaSBarry Smith   For two tuples x and y with the same degree sum, partial degree sums over the final elements of the tuples
957fbdc3dfeSToby Isaac   acts as a tiebreaker.  For example, (2, 1, 1) and (1, 2, 1) have the same degree sum, but the degree sum over the
958fbdc3dfeSToby Isaac   last two elements is smaller for the former, so (2, 1, 1) < (1, 2, 1).
959fbdc3dfeSToby Isaac 
960db781477SPatrick Sanan .seealso: `PetscDTIndexToGradedOrder()`
961fbdc3dfeSToby Isaac @*/
962d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTGradedOrderToIndex(PetscInt len, const PetscInt degtup[], PetscInt *index)
963d71ae5a4SJacob Faibussowitsch {
964fbdc3dfeSToby Isaac   PetscInt i, idx, sum, total;
965fbdc3dfeSToby Isaac 
966fbdc3dfeSToby Isaac   PetscFunctionBeginHot;
96708401ef6SPierre Jolivet   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
968fbdc3dfeSToby Isaac   for (i = 0, sum = 0; i < len; i++) sum += degtup[i];
969fbdc3dfeSToby Isaac   idx   = 0;
970fbdc3dfeSToby Isaac   total = 1;
971fbdc3dfeSToby Isaac   for (i = 0; i < sum; i++) {
972fbdc3dfeSToby Isaac     idx += total;
973fbdc3dfeSToby Isaac     total = (total * (len + i)) / (i + 1);
974fbdc3dfeSToby Isaac   }
975fbdc3dfeSToby Isaac   for (i = 0; i < len - 1; i++) {
976fbdc3dfeSToby Isaac     PetscInt c;
977fbdc3dfeSToby Isaac 
978fbdc3dfeSToby Isaac     total = 1;
979fbdc3dfeSToby Isaac     sum -= degtup[i];
980fbdc3dfeSToby Isaac     for (c = 0; c < sum; c++) {
981fbdc3dfeSToby Isaac       idx += total;
982fbdc3dfeSToby Isaac       total = (total * (len - 1 - i + c)) / (c + 1);
983fbdc3dfeSToby Isaac     }
984fbdc3dfeSToby Isaac   }
985fbdc3dfeSToby Isaac   *index = idx;
9863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
987fbdc3dfeSToby Isaac }
988fbdc3dfeSToby Isaac 
989e3aa2e09SToby Isaac static PetscBool PKDCite       = PETSC_FALSE;
990e3aa2e09SToby Isaac const char       PKDCitation[] = "@article{Kirby2010,\n"
991e3aa2e09SToby Isaac                                  "  title={Singularity-free evaluation of collapsed-coordinate orthogonal polynomials},\n"
992e3aa2e09SToby Isaac                                  "  author={Kirby, Robert C},\n"
993e3aa2e09SToby Isaac                                  "  journal={ACM Transactions on Mathematical Software (TOMS)},\n"
994e3aa2e09SToby Isaac                                  "  volume={37},\n"
995e3aa2e09SToby Isaac                                  "  number={1},\n"
996e3aa2e09SToby Isaac                                  "  pages={1--16},\n"
997e3aa2e09SToby Isaac                                  "  year={2010},\n"
998e3aa2e09SToby Isaac                                  "  publisher={ACM New York, NY, USA}\n}\n";
999e3aa2e09SToby Isaac 
1000fbdc3dfeSToby Isaac /*@
1001d8f25ad8SToby Isaac   PetscDTPKDEvalJet - Evaluate the jet (function and derivatives) of the Proriol-Koornwinder-Dubiner (PKD) basis for
1002fbdc3dfeSToby Isaac   the space of polynomials up to a given degree.  The PKD basis is L2-orthonormal on the biunit simplex (which is used
1003fbdc3dfeSToby Isaac   as the reference element for finite elements in PETSc), which makes it a stable basis to use for evaluating
1004fbdc3dfeSToby Isaac   polynomials in that domain.
1005fbdc3dfeSToby Isaac 
10064165533cSJose E. Roman   Input Parameters:
1007fbdc3dfeSToby Isaac + dim - the number of variables in the multivariate polynomials
1008fbdc3dfeSToby Isaac . npoints - the number of points to evaluate the polynomials at
1009fbdc3dfeSToby Isaac . points - [npoints x dim] array of point coordinates
1010fbdc3dfeSToby Isaac . 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.
1011fbdc3dfeSToby Isaac - k - the maximum order partial derivative to evaluate in the jet.  There are (dim + k choose dim) partial derivatives
1012fbdc3dfeSToby Isaac   in the jet.  Choosing k = 0 means to evaluate just the function and no derivatives
1013fbdc3dfeSToby Isaac 
10146aad120cSJose E. Roman   Output Parameters:
1015fbdc3dfeSToby Isaac - p - an array containing the evaluations of the PKD polynomials' jets on the points.  The size is ((dim + degree)
1016fbdc3dfeSToby Isaac   choose dim) x ((dim + k) choose dim) x npoints, which also describes the order of the dimensions of this
1017fbdc3dfeSToby Isaac   three-dimensional array: the first (slowest varying) dimension is basis function index; the second dimension is jet
1018fbdc3dfeSToby Isaac   index; the third (fastest varying) dimension is the index of the evaluation point.
1019fbdc3dfeSToby Isaac 
1020fbdc3dfeSToby Isaac   Level: advanced
1021fbdc3dfeSToby Isaac 
1022dce8aebaSBarry Smith   Notes:
1023dce8aebaSBarry Smith   The ordering of the basis functions, and the ordering of the derivatives in the jet, both follow the graded
1024dce8aebaSBarry Smith   ordering of `PetscDTIndexToGradedOrder()` and `PetscDTGradedOrderToIndex()`.  For example, in 3D, the polynomial with
1025dce8aebaSBarry Smith   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);
1026fbdc3dfeSToby Isaac   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).
1027fbdc3dfeSToby Isaac 
1028e3aa2e09SToby Isaac   The implementation uses Kirby's singularity-free evaluation algorithm, https://doi.org/10.1145/1644001.1644006.
1029e3aa2e09SToby Isaac 
1030db781477SPatrick Sanan .seealso: `PetscDTGradedOrderToIndex()`, `PetscDTIndexToGradedOrder()`, `PetscDTJacobiEvalJet()`
1031fbdc3dfeSToby Isaac @*/
1032d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTPKDEvalJet(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt k, PetscReal p[])
1033d71ae5a4SJacob Faibussowitsch {
1034fbdc3dfeSToby Isaac   PetscInt   degidx, kidx, d, pt;
1035fbdc3dfeSToby Isaac   PetscInt   Nk, Ndeg;
1036fbdc3dfeSToby Isaac   PetscInt  *ktup, *degtup;
1037fbdc3dfeSToby Isaac   PetscReal *scales, initscale, scaleexp;
1038fbdc3dfeSToby Isaac 
1039fbdc3dfeSToby Isaac   PetscFunctionBegin;
10409566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(PKDCitation, &PKDCite));
10419566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + k, k, &Nk));
10429566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(degree + dim, degree, &Ndeg));
10439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim, &degtup, dim, &ktup));
10449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Ndeg, &scales));
1045fbdc3dfeSToby Isaac   initscale = 1.;
1046fbdc3dfeSToby Isaac   if (dim > 1) {
10479566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomial(dim, 2, &scaleexp));
10482f613bf5SBarry Smith     initscale = PetscPowReal(2., scaleexp * 0.5);
1049fbdc3dfeSToby Isaac   }
1050fbdc3dfeSToby Isaac   for (degidx = 0; degidx < Ndeg; degidx++) {
1051fbdc3dfeSToby Isaac     PetscInt  e, i;
1052fbdc3dfeSToby Isaac     PetscInt  m1idx = -1, m2idx = -1;
1053fbdc3dfeSToby Isaac     PetscInt  n;
1054fbdc3dfeSToby Isaac     PetscInt  degsum;
1055fbdc3dfeSToby Isaac     PetscReal alpha;
1056fbdc3dfeSToby Isaac     PetscReal cnm1, cnm1x, cnm2;
1057fbdc3dfeSToby Isaac     PetscReal norm;
1058fbdc3dfeSToby Isaac 
10599566063dSJacob Faibussowitsch     PetscCall(PetscDTIndexToGradedOrder(dim, degidx, degtup));
10609371c9d4SSatish Balay     for (d = dim - 1; d >= 0; d--)
10619371c9d4SSatish Balay       if (degtup[d]) break;
1062fbdc3dfeSToby Isaac     if (d < 0) { /* constant is 1 everywhere, all derivatives are zero */
1063fbdc3dfeSToby Isaac       scales[degidx] = initscale;
1064fbdc3dfeSToby Isaac       for (e = 0; e < dim; e++) {
10659566063dSJacob Faibussowitsch         PetscCall(PetscDTJacobiNorm(e, 0., 0, &norm));
1066fbdc3dfeSToby Isaac         scales[degidx] /= norm;
1067fbdc3dfeSToby Isaac       }
1068fbdc3dfeSToby Isaac       for (i = 0; i < npoints; i++) p[degidx * Nk * npoints + i] = 1.;
1069fbdc3dfeSToby Isaac       for (i = 0; i < (Nk - 1) * npoints; i++) p[(degidx * Nk + 1) * npoints + i] = 0.;
1070fbdc3dfeSToby Isaac       continue;
1071fbdc3dfeSToby Isaac     }
1072fbdc3dfeSToby Isaac     n = degtup[d];
1073fbdc3dfeSToby Isaac     degtup[d]--;
10749566063dSJacob Faibussowitsch     PetscCall(PetscDTGradedOrderToIndex(dim, degtup, &m1idx));
1075fbdc3dfeSToby Isaac     if (degtup[d] > 0) {
1076fbdc3dfeSToby Isaac       degtup[d]--;
10779566063dSJacob Faibussowitsch       PetscCall(PetscDTGradedOrderToIndex(dim, degtup, &m2idx));
1078fbdc3dfeSToby Isaac       degtup[d]++;
1079fbdc3dfeSToby Isaac     }
1080fbdc3dfeSToby Isaac     degtup[d]++;
1081fbdc3dfeSToby Isaac     for (e = 0, degsum = 0; e < d; e++) degsum += degtup[e];
1082fbdc3dfeSToby Isaac     alpha = 2 * degsum + d;
1083fbdc3dfeSToby Isaac     PetscDTJacobiRecurrence_Internal(n, alpha, 0., cnm1, cnm1x, cnm2);
1084fbdc3dfeSToby Isaac 
1085fbdc3dfeSToby Isaac     scales[degidx] = initscale;
1086fbdc3dfeSToby Isaac     for (e = 0, degsum = 0; e < dim; e++) {
1087fbdc3dfeSToby Isaac       PetscInt  f;
1088fbdc3dfeSToby Isaac       PetscReal ealpha;
1089fbdc3dfeSToby Isaac       PetscReal enorm;
1090fbdc3dfeSToby Isaac 
1091fbdc3dfeSToby Isaac       ealpha = 2 * degsum + e;
1092fbdc3dfeSToby Isaac       for (f = 0; f < degsum; f++) scales[degidx] *= 2.;
10939566063dSJacob Faibussowitsch       PetscCall(PetscDTJacobiNorm(ealpha, 0., degtup[e], &enorm));
1094fbdc3dfeSToby Isaac       scales[degidx] /= enorm;
1095fbdc3dfeSToby Isaac       degsum += degtup[e];
1096fbdc3dfeSToby Isaac     }
1097fbdc3dfeSToby Isaac 
1098fbdc3dfeSToby Isaac     for (pt = 0; pt < npoints; pt++) {
1099fbdc3dfeSToby Isaac       /* compute the multipliers */
1100fbdc3dfeSToby Isaac       PetscReal thetanm1, thetanm1x, thetanm2;
1101fbdc3dfeSToby Isaac 
1102fbdc3dfeSToby Isaac       thetanm1x = dim - (d + 1) + 2. * points[pt * dim + d];
1103fbdc3dfeSToby Isaac       for (e = d + 1; e < dim; e++) thetanm1x += points[pt * dim + e];
1104fbdc3dfeSToby Isaac       thetanm1x *= 0.5;
1105fbdc3dfeSToby Isaac       thetanm1 = (2. - (dim - (d + 1)));
1106fbdc3dfeSToby Isaac       for (e = d + 1; e < dim; e++) thetanm1 -= points[pt * dim + e];
1107fbdc3dfeSToby Isaac       thetanm1 *= 0.5;
1108fbdc3dfeSToby Isaac       thetanm2 = thetanm1 * thetanm1;
1109fbdc3dfeSToby Isaac 
1110fbdc3dfeSToby Isaac       for (kidx = 0; kidx < Nk; kidx++) {
1111fbdc3dfeSToby Isaac         PetscInt f;
1112fbdc3dfeSToby Isaac 
11139566063dSJacob Faibussowitsch         PetscCall(PetscDTIndexToGradedOrder(dim, kidx, ktup));
1114fbdc3dfeSToby Isaac         /* first sum in the same derivative terms */
1115fbdc3dfeSToby Isaac         p[(degidx * Nk + kidx) * npoints + pt] = (cnm1 * thetanm1 + cnm1x * thetanm1x) * p[(m1idx * Nk + kidx) * npoints + pt];
1116ad540459SPierre Jolivet         if (m2idx >= 0) p[(degidx * Nk + kidx) * npoints + pt] -= cnm2 * thetanm2 * p[(m2idx * Nk + kidx) * npoints + pt];
1117fbdc3dfeSToby Isaac 
1118fbdc3dfeSToby Isaac         for (f = d; f < dim; f++) {
1119fbdc3dfeSToby Isaac           PetscInt km1idx, mplty = ktup[f];
1120fbdc3dfeSToby Isaac 
1121fbdc3dfeSToby Isaac           if (!mplty) continue;
1122fbdc3dfeSToby Isaac           ktup[f]--;
11239566063dSJacob Faibussowitsch           PetscCall(PetscDTGradedOrderToIndex(dim, ktup, &km1idx));
1124fbdc3dfeSToby Isaac 
1125fbdc3dfeSToby Isaac           /* the derivative of  cnm1x * thetanm1x  wrt x variable f is 0.5 * cnm1x if f > d otherwise it is cnm1x */
1126fbdc3dfeSToby Isaac           /* the derivative of  cnm1  * thetanm1   wrt x variable f is 0 if f == d, otherwise it is -0.5 * cnm1 */
1127fbdc3dfeSToby Isaac           /* the derivative of -cnm2  * thetanm2   wrt x variable f is 0 if f == d, otherwise it is cnm2 * thetanm1 */
1128fbdc3dfeSToby Isaac           if (f > d) {
1129fbdc3dfeSToby Isaac             PetscInt f2;
1130fbdc3dfeSToby Isaac 
1131fbdc3dfeSToby Isaac             p[(degidx * Nk + kidx) * npoints + pt] += mplty * 0.5 * (cnm1x - cnm1) * p[(m1idx * Nk + km1idx) * npoints + pt];
1132fbdc3dfeSToby Isaac             if (m2idx >= 0) {
1133fbdc3dfeSToby Isaac               p[(degidx * Nk + kidx) * npoints + pt] += mplty * cnm2 * thetanm1 * p[(m2idx * Nk + km1idx) * npoints + pt];
1134fbdc3dfeSToby Isaac               /* second derivatives of -cnm2  * thetanm2   wrt x variable f,f2 is like - 0.5 * cnm2 */
1135fbdc3dfeSToby Isaac               for (f2 = f; f2 < dim; f2++) {
1136fbdc3dfeSToby Isaac                 PetscInt km2idx, mplty2 = ktup[f2];
1137fbdc3dfeSToby Isaac                 PetscInt factor;
1138fbdc3dfeSToby Isaac 
1139fbdc3dfeSToby Isaac                 if (!mplty2) continue;
1140fbdc3dfeSToby Isaac                 ktup[f2]--;
11419566063dSJacob Faibussowitsch                 PetscCall(PetscDTGradedOrderToIndex(dim, ktup, &km2idx));
1142fbdc3dfeSToby Isaac 
1143fbdc3dfeSToby Isaac                 factor = mplty * mplty2;
1144fbdc3dfeSToby Isaac                 if (f == f2) factor /= 2;
1145fbdc3dfeSToby Isaac                 p[(degidx * Nk + kidx) * npoints + pt] -= 0.5 * factor * cnm2 * p[(m2idx * Nk + km2idx) * npoints + pt];
1146fbdc3dfeSToby Isaac                 ktup[f2]++;
1147fbdc3dfeSToby Isaac               }
11483034baaeSToby Isaac             }
1149fbdc3dfeSToby Isaac           } else {
1150fbdc3dfeSToby Isaac             p[(degidx * Nk + kidx) * npoints + pt] += mplty * cnm1x * p[(m1idx * Nk + km1idx) * npoints + pt];
1151fbdc3dfeSToby Isaac           }
1152fbdc3dfeSToby Isaac           ktup[f]++;
1153fbdc3dfeSToby Isaac         }
1154fbdc3dfeSToby Isaac       }
1155fbdc3dfeSToby Isaac     }
1156fbdc3dfeSToby Isaac   }
1157fbdc3dfeSToby Isaac   for (degidx = 0; degidx < Ndeg; degidx++) {
1158fbdc3dfeSToby Isaac     PetscReal scale = scales[degidx];
1159fbdc3dfeSToby Isaac     PetscInt  i;
1160fbdc3dfeSToby Isaac 
1161fbdc3dfeSToby Isaac     for (i = 0; i < Nk * npoints; i++) p[degidx * Nk * npoints + i] *= scale;
1162fbdc3dfeSToby Isaac   }
11639566063dSJacob Faibussowitsch   PetscCall(PetscFree(scales));
11649566063dSJacob Faibussowitsch   PetscCall(PetscFree2(degtup, ktup));
11653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1166fbdc3dfeSToby Isaac }
1167fbdc3dfeSToby Isaac 
1168d8f25ad8SToby Isaac /*@
1169d8f25ad8SToby Isaac   PetscDTPTrimmedSize - The size of the trimmed polynomial space of k-forms with a given degree and form degree,
1170dce8aebaSBarry Smith   which can be evaluated in `PetscDTPTrimmedEvalJet()`.
1171d8f25ad8SToby Isaac 
1172d8f25ad8SToby Isaac   Input Parameters:
1173d8f25ad8SToby Isaac + dim - the number of variables in the multivariate polynomials
1174d8f25ad8SToby Isaac . degree - the degree (sum of degrees on the variables in a monomial) of the trimmed polynomial space.
1175d8f25ad8SToby Isaac - formDegree - the degree of the form
1176d8f25ad8SToby Isaac 
11776aad120cSJose E. Roman   Output Parameters:
117820f4b53cSBarry Smith - size - The number ((`dim` + `degree`) choose (`dim` + `formDegree`)) x ((`degree` + `formDegree` - 1) choose (`formDegree`))
1179d8f25ad8SToby Isaac 
1180d8f25ad8SToby Isaac   Level: advanced
1181d8f25ad8SToby Isaac 
1182db781477SPatrick Sanan .seealso: `PetscDTPTrimmedEvalJet()`
1183d8f25ad8SToby Isaac @*/
1184d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTPTrimmedSize(PetscInt dim, PetscInt degree, PetscInt formDegree, PetscInt *size)
1185d71ae5a4SJacob Faibussowitsch {
1186d8f25ad8SToby Isaac   PetscInt Nrk, Nbpt; // number of trimmed polynomials
1187d8f25ad8SToby Isaac 
1188d8f25ad8SToby Isaac   PetscFunctionBegin;
1189d8f25ad8SToby Isaac   formDegree = PetscAbsInt(formDegree);
11909566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(degree + dim, degree + formDegree, &Nbpt));
11919566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(degree + formDegree - 1, formDegree, &Nrk));
1192d8f25ad8SToby Isaac   Nbpt *= Nrk;
1193d8f25ad8SToby Isaac   *size = Nbpt;
11943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1195d8f25ad8SToby Isaac }
1196d8f25ad8SToby Isaac 
1197d8f25ad8SToby Isaac /* there was a reference implementation based on section 4.4 of Arnold, Falk & Winther (acta numerica, 2006), but it
1198d8f25ad8SToby Isaac  * was inferior to this implementation */
1199d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTPTrimmedEvalJet_Internal(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt formDegree, PetscInt jetDegree, PetscReal p[])
1200d71ae5a4SJacob Faibussowitsch {
1201d8f25ad8SToby Isaac   PetscInt  formDegreeOrig = formDegree;
1202d8f25ad8SToby Isaac   PetscBool formNegative   = (formDegreeOrig < 0) ? PETSC_TRUE : PETSC_FALSE;
1203d8f25ad8SToby Isaac 
1204d8f25ad8SToby Isaac   PetscFunctionBegin;
1205d8f25ad8SToby Isaac   formDegree = PetscAbsInt(formDegreeOrig);
1206d8f25ad8SToby Isaac   if (formDegree == 0) {
12079566063dSJacob Faibussowitsch     PetscCall(PetscDTPKDEvalJet(dim, npoints, points, degree, jetDegree, p));
12083ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1209d8f25ad8SToby Isaac   }
1210d8f25ad8SToby Isaac   if (formDegree == dim) {
12119566063dSJacob Faibussowitsch     PetscCall(PetscDTPKDEvalJet(dim, npoints, points, degree - 1, jetDegree, p));
12123ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1213d8f25ad8SToby Isaac   }
1214d8f25ad8SToby Isaac   PetscInt Nbpt;
12159566063dSJacob Faibussowitsch   PetscCall(PetscDTPTrimmedSize(dim, degree, formDegree, &Nbpt));
1216d8f25ad8SToby Isaac   PetscInt Nf;
12179566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, formDegree, &Nf));
1218d8f25ad8SToby Isaac   PetscInt Nk;
12199566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + jetDegree, dim, &Nk));
12209566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(p, Nbpt * Nf * Nk * npoints));
1221d8f25ad8SToby Isaac 
1222d8f25ad8SToby Isaac   PetscInt Nbpm1; // number of scalar polynomials up to degree - 1;
12239566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + degree - 1, dim, &Nbpm1));
1224d8f25ad8SToby Isaac   PetscReal *p_scalar;
12259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nbpm1 * Nk * npoints, &p_scalar));
12269566063dSJacob Faibussowitsch   PetscCall(PetscDTPKDEvalJet(dim, npoints, points, degree - 1, jetDegree, p_scalar));
1227d8f25ad8SToby Isaac   PetscInt total = 0;
1228d8f25ad8SToby Isaac   // First add the full polynomials up to degree - 1 into the basis: take the scalar
1229d8f25ad8SToby Isaac   // and copy one for each form component
1230d8f25ad8SToby Isaac   for (PetscInt i = 0; i < Nbpm1; i++) {
1231d8f25ad8SToby Isaac     const PetscReal *src = &p_scalar[i * Nk * npoints];
1232d8f25ad8SToby Isaac     for (PetscInt f = 0; f < Nf; f++) {
1233d8f25ad8SToby Isaac       PetscReal *dest = &p[(total++ * Nf + f) * Nk * npoints];
12349566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(dest, src, Nk * npoints));
1235d8f25ad8SToby Isaac     }
1236d8f25ad8SToby Isaac   }
1237d8f25ad8SToby Isaac   PetscInt *form_atoms;
12389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(formDegree + 1, &form_atoms));
1239d8f25ad8SToby Isaac   // construct the interior product pattern
1240d8f25ad8SToby Isaac   PetscInt(*pattern)[3];
1241d8f25ad8SToby Isaac   PetscInt Nf1; // number of formDegree + 1 forms
12429566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, formDegree + 1, &Nf1));
1243d8f25ad8SToby Isaac   PetscInt nnz = Nf1 * (formDegree + 1);
12449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf1 * (formDegree + 1), &pattern));
12459566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVInteriorPattern(dim, formDegree + 1, pattern));
1246d8f25ad8SToby Isaac   PetscReal centroid = (1. - dim) / (dim + 1.);
1247d8f25ad8SToby Isaac   PetscInt *deriv;
12489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim, &deriv));
1249d8f25ad8SToby Isaac   for (PetscInt d = dim; d >= formDegree + 1; d--) {
1250d8f25ad8SToby Isaac     PetscInt Nfd1; // number of formDegree + 1 forms in dimension d that include dx_0
1251d8f25ad8SToby Isaac                    // (equal to the number of formDegree forms in dimension d-1)
12529566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(d - 1, formDegree, &Nfd1));
1253d8f25ad8SToby Isaac     // The number of homogeneous (degree-1) scalar polynomials in d variables
1254d8f25ad8SToby Isaac     PetscInt Nh;
12559566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(d - 1 + degree - 1, d - 1, &Nh));
1256d8f25ad8SToby Isaac     const PetscReal *h_scalar = &p_scalar[(Nbpm1 - Nh) * Nk * npoints];
1257d8f25ad8SToby Isaac     for (PetscInt b = 0; b < Nh; b++) {
1258d8f25ad8SToby Isaac       const PetscReal *h_s = &h_scalar[b * Nk * npoints];
1259d8f25ad8SToby Isaac       for (PetscInt f = 0; f < Nfd1; f++) {
1260d8f25ad8SToby Isaac         // construct all formDegree+1 forms that start with dx_(dim - d) /\ ...
1261d8f25ad8SToby Isaac         form_atoms[0] = dim - d;
12629566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSubset(d - 1, formDegree, f, &form_atoms[1]));
1263ad540459SPierre Jolivet         for (PetscInt i = 0; i < formDegree; i++) form_atoms[1 + i] += form_atoms[0] + 1;
1264d8f25ad8SToby Isaac         PetscInt f_ind; // index of the resulting form
12659566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(dim, formDegree + 1, form_atoms, &f_ind));
1266d8f25ad8SToby Isaac         PetscReal *p_f = &p[total++ * Nf * Nk * npoints];
1267d8f25ad8SToby Isaac         for (PetscInt nz = 0; nz < nnz; nz++) {
1268d8f25ad8SToby Isaac           PetscInt  i     = pattern[nz][0]; // formDegree component
1269d8f25ad8SToby Isaac           PetscInt  j     = pattern[nz][1]; // (formDegree + 1) component
1270d8f25ad8SToby Isaac           PetscInt  v     = pattern[nz][2]; // coordinate component
1271d8f25ad8SToby Isaac           PetscReal scale = v < 0 ? -1. : 1.;
1272d8f25ad8SToby Isaac 
1273d8f25ad8SToby Isaac           i     = formNegative ? (Nf - 1 - i) : i;
1274d8f25ad8SToby Isaac           scale = (formNegative && (i & 1)) ? -scale : scale;
1275d8f25ad8SToby Isaac           v     = v < 0 ? -(v + 1) : v;
1276ad540459SPierre Jolivet           if (j != f_ind) continue;
1277d8f25ad8SToby Isaac           PetscReal *p_i = &p_f[i * Nk * npoints];
1278d8f25ad8SToby Isaac           for (PetscInt jet = 0; jet < Nk; jet++) {
1279d8f25ad8SToby Isaac             const PetscReal *h_jet = &h_s[jet * npoints];
1280d8f25ad8SToby Isaac             PetscReal       *p_jet = &p_i[jet * npoints];
1281d8f25ad8SToby Isaac 
1282ad540459SPierre Jolivet             for (PetscInt pt = 0; pt < npoints; pt++) p_jet[pt] += scale * h_jet[pt] * (points[pt * dim + v] - centroid);
12839566063dSJacob Faibussowitsch             PetscCall(PetscDTIndexToGradedOrder(dim, jet, deriv));
1284d8f25ad8SToby Isaac             deriv[v]++;
1285d8f25ad8SToby Isaac             PetscReal mult = deriv[v];
1286d8f25ad8SToby Isaac             PetscInt  l;
12879566063dSJacob Faibussowitsch             PetscCall(PetscDTGradedOrderToIndex(dim, deriv, &l));
1288ad540459SPierre Jolivet             if (l >= Nk) continue;
1289d8f25ad8SToby Isaac             p_jet = &p_i[l * npoints];
1290ad540459SPierre Jolivet             for (PetscInt pt = 0; pt < npoints; pt++) p_jet[pt] += scale * mult * h_jet[pt];
1291d8f25ad8SToby Isaac             deriv[v]--;
1292d8f25ad8SToby Isaac           }
1293d8f25ad8SToby Isaac         }
1294d8f25ad8SToby Isaac       }
1295d8f25ad8SToby Isaac     }
1296d8f25ad8SToby Isaac   }
129708401ef6SPierre Jolivet   PetscCheck(total == Nbpt, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrectly counted P trimmed polynomials");
12989566063dSJacob Faibussowitsch   PetscCall(PetscFree(deriv));
12999566063dSJacob Faibussowitsch   PetscCall(PetscFree(pattern));
13009566063dSJacob Faibussowitsch   PetscCall(PetscFree(form_atoms));
13019566063dSJacob Faibussowitsch   PetscCall(PetscFree(p_scalar));
13023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1303d8f25ad8SToby Isaac }
1304d8f25ad8SToby Isaac 
1305d8f25ad8SToby Isaac /*@
1306d8f25ad8SToby Isaac   PetscDTPTrimmedEvalJet - Evaluate the jet (function and derivatives) of a basis of the trimmed polynomial k-forms up to
1307d8f25ad8SToby Isaac   a given degree.
1308d8f25ad8SToby Isaac 
1309d8f25ad8SToby Isaac   Input Parameters:
1310d8f25ad8SToby Isaac + dim - the number of variables in the multivariate polynomials
1311d8f25ad8SToby Isaac . npoints - the number of points to evaluate the polynomials at
1312d8f25ad8SToby Isaac . points - [npoints x dim] array of point coordinates
1313d8f25ad8SToby Isaac . degree - the degree (sum of degrees on the variables in a monomial) of the trimmed polynomial space to evaluate.
1314d8f25ad8SToby Isaac            There are ((dim + degree) choose (dim + formDegree)) x ((degree + formDegree - 1) choose (formDegree)) polynomials in this space.
1315dce8aebaSBarry Smith            (You can use `PetscDTPTrimmedSize()` to compute this size.)
1316d8f25ad8SToby Isaac . formDegree - the degree of the form
1317d8f25ad8SToby Isaac - jetDegree - the maximum order partial derivative to evaluate in the jet.  There are ((dim + jetDegree) choose dim) partial derivatives
1318d8f25ad8SToby Isaac               in the jet.  Choosing jetDegree = 0 means to evaluate just the function and no derivatives
1319d8f25ad8SToby Isaac 
132020f4b53cSBarry Smith   Output Parameter:
132120f4b53cSBarry Smith . p - an array containing the evaluations of the PKD polynomials' jets on the points.  The size is
1322dce8aebaSBarry Smith       `PetscDTPTrimmedSize()` x ((dim + formDegree) choose dim) x ((dim + k) choose dim) x npoints,
1323d8f25ad8SToby Isaac       which also describes the order of the dimensions of this
1324d8f25ad8SToby Isaac       four-dimensional array:
1325d8f25ad8SToby Isaac         the first (slowest varying) dimension is basis function index;
1326d8f25ad8SToby Isaac         the second dimension is component of the form;
1327d8f25ad8SToby Isaac         the third dimension is jet index;
1328d8f25ad8SToby Isaac         the fourth (fastest varying) dimension is the index of the evaluation point.
1329d8f25ad8SToby Isaac 
1330d8f25ad8SToby Isaac   Level: advanced
1331d8f25ad8SToby Isaac 
1332dce8aebaSBarry Smith   Notes:
1333dce8aebaSBarry Smith   The ordering of the basis functions is not graded, so the basis functions are not nested by degree like `PetscDTPKDEvalJet()`.
1334d8f25ad8SToby Isaac   The basis functions are not an L2-orthonormal basis on any particular domain.
1335d8f25ad8SToby Isaac 
1336d8f25ad8SToby Isaac   The implementation is based on the description of the trimmed polynomials up to degree r as
1337d8f25ad8SToby Isaac   the direct sum of polynomials up to degree (r-1) and the Koszul differential applied to
1338d8f25ad8SToby Isaac   homogeneous polynomials of degree (r-1).
1339d8f25ad8SToby Isaac 
1340db781477SPatrick Sanan .seealso: `PetscDTPKDEvalJet()`, `PetscDTPTrimmedSize()`
1341d8f25ad8SToby Isaac @*/
1342d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt formDegree, PetscInt jetDegree, PetscReal p[])
1343d71ae5a4SJacob Faibussowitsch {
1344d8f25ad8SToby Isaac   PetscFunctionBegin;
13459566063dSJacob Faibussowitsch   PetscCall(PetscDTPTrimmedEvalJet_Internal(dim, npoints, points, degree, formDegree, jetDegree, p));
13463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1347d8f25ad8SToby Isaac }
1348d8f25ad8SToby Isaac 
1349e6a796c3SToby Isaac /* solve the symmetric tridiagonal eigenvalue system, writing the eigenvalues into eigs and the eigenvectors into V
1350e6a796c3SToby Isaac  * with lds n; diag and subdiag are overwritten */
1351d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTSymmetricTridiagonalEigensolve(PetscInt n, PetscReal diag[], PetscReal subdiag[], PetscReal eigs[], PetscScalar V[])
1352d71ae5a4SJacob Faibussowitsch {
1353e6a796c3SToby Isaac   char          jobz   = 'V'; /* eigenvalues and eigenvectors */
1354e6a796c3SToby Isaac   char          range  = 'A'; /* all eigenvalues will be found */
1355e6a796c3SToby Isaac   PetscReal     VL     = 0.;  /* ignored because range is 'A' */
1356e6a796c3SToby Isaac   PetscReal     VU     = 0.;  /* ignored because range is 'A' */
1357e6a796c3SToby Isaac   PetscBLASInt  IL     = 0;   /* ignored because range is 'A' */
1358e6a796c3SToby Isaac   PetscBLASInt  IU     = 0;   /* ignored because range is 'A' */
1359e6a796c3SToby Isaac   PetscReal     abstol = 0.;  /* unused */
1360e6a796c3SToby Isaac   PetscBLASInt  bn, bm, ldz;  /* bm will equal bn on exit */
1361e6a796c3SToby Isaac   PetscBLASInt *isuppz;
1362e6a796c3SToby Isaac   PetscBLASInt  lwork, liwork;
1363e6a796c3SToby Isaac   PetscReal     workquery;
1364e6a796c3SToby Isaac   PetscBLASInt  iworkquery;
1365e6a796c3SToby Isaac   PetscBLASInt *iwork;
1366e6a796c3SToby Isaac   PetscBLASInt  info;
1367e6a796c3SToby Isaac   PetscReal    *work = NULL;
1368e6a796c3SToby Isaac 
1369e6a796c3SToby Isaac   PetscFunctionBegin;
1370e6a796c3SToby Isaac #if !defined(PETSCDTGAUSSIANQUADRATURE_EIG)
1371e6a796c3SToby Isaac   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
1372e6a796c3SToby Isaac #endif
13739566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &bn));
13749566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &ldz));
1375e6a796c3SToby Isaac #if !defined(PETSC_MISSING_LAPACK_STEGR)
13769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * n, &isuppz));
1377e6a796c3SToby Isaac   lwork  = -1;
1378e6a796c3SToby Isaac   liwork = -1;
1379792fecdfSBarry Smith   PetscCallBLAS("LAPACKstegr", LAPACKstegr_(&jobz, &range, &bn, diag, subdiag, &VL, &VU, &IL, &IU, &abstol, &bm, eigs, V, &ldz, isuppz, &workquery, &lwork, &iworkquery, &liwork, &info));
138028b400f6SJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xSTEGR error");
1381e6a796c3SToby Isaac   lwork  = (PetscBLASInt)workquery;
1382e6a796c3SToby Isaac   liwork = (PetscBLASInt)iworkquery;
13839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(lwork, &work, liwork, &iwork));
13849566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1385792fecdfSBarry Smith   PetscCallBLAS("LAPACKstegr", LAPACKstegr_(&jobz, &range, &bn, diag, subdiag, &VL, &VU, &IL, &IU, &abstol, &bm, eigs, V, &ldz, isuppz, work, &lwork, iwork, &liwork, &info));
13869566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPop());
138728b400f6SJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xSTEGR error");
13889566063dSJacob Faibussowitsch   PetscCall(PetscFree2(work, iwork));
13899566063dSJacob Faibussowitsch   PetscCall(PetscFree(isuppz));
1390e6a796c3SToby Isaac #elif !defined(PETSC_MISSING_LAPACK_STEQR)
1391e6a796c3SToby Isaac   jobz = 'I'; /* Compute eigenvalues and eigenvectors of the
1392e6a796c3SToby Isaac                  tridiagonal matrix.  Z is initialized to the identity
1393e6a796c3SToby Isaac                  matrix. */
13949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(PetscMax(1, 2 * n - 2), &work));
1395792fecdfSBarry Smith   PetscCallBLAS("LAPACKsteqr", LAPACKsteqr_("I", &bn, diag, subdiag, V, &ldz, work, &info));
13969566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPop());
139728b400f6SJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xSTEQR error");
13989566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
13999566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(eigs, diag, n));
1400e6a796c3SToby Isaac #endif
14013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1402e6a796c3SToby Isaac }
1403e6a796c3SToby Isaac 
1404e6a796c3SToby Isaac /* Formula for the weights at the endpoints (-1 and 1) of Gauss-Lobatto-Jacobi
1405e6a796c3SToby Isaac  * quadrature rules on the interval [-1, 1] */
1406d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTGaussLobattoJacobiEndweights_Internal(PetscInt n, PetscReal alpha, PetscReal beta, PetscReal *leftw, PetscReal *rightw)
1407d71ae5a4SJacob Faibussowitsch {
1408e6a796c3SToby Isaac   PetscReal twoab1;
1409e6a796c3SToby Isaac   PetscInt  m = n - 2;
1410e6a796c3SToby Isaac   PetscReal a = alpha + 1.;
1411e6a796c3SToby Isaac   PetscReal b = beta + 1.;
1412e6a796c3SToby Isaac   PetscReal gra, grb;
1413e6a796c3SToby Isaac 
1414e6a796c3SToby Isaac   PetscFunctionBegin;
1415e6a796c3SToby Isaac   twoab1 = PetscPowReal(2., a + b - 1.);
1416e6a796c3SToby Isaac #if defined(PETSC_HAVE_LGAMMA)
14179371c9d4SSatish Balay   grb = PetscExpReal(2. * PetscLGamma(b + 1.) + PetscLGamma(m + 1.) + PetscLGamma(m + a + 1.) - (PetscLGamma(m + b + 1) + PetscLGamma(m + a + b + 1.)));
14189371c9d4SSatish Balay   gra = PetscExpReal(2. * PetscLGamma(a + 1.) + PetscLGamma(m + 1.) + PetscLGamma(m + b + 1.) - (PetscLGamma(m + a + 1) + PetscLGamma(m + a + b + 1.)));
1419e6a796c3SToby Isaac #else
1420e6a796c3SToby Isaac   {
1421e6a796c3SToby Isaac     PetscInt alphai = (PetscInt)alpha;
1422e6a796c3SToby Isaac     PetscInt betai  = (PetscInt)beta;
1423e6a796c3SToby Isaac 
1424e6a796c3SToby Isaac     if ((PetscReal)alphai == alpha && (PetscReal)betai == beta) {
1425e6a796c3SToby Isaac       PetscReal binom1, binom2;
1426e6a796c3SToby Isaac 
14279566063dSJacob Faibussowitsch       PetscCall(PetscDTBinomial(m + b, b, &binom1));
14289566063dSJacob Faibussowitsch       PetscCall(PetscDTBinomial(m + a + b, b, &binom2));
1429e6a796c3SToby Isaac       grb = 1. / (binom1 * binom2);
14309566063dSJacob Faibussowitsch       PetscCall(PetscDTBinomial(m + a, a, &binom1));
14319566063dSJacob Faibussowitsch       PetscCall(PetscDTBinomial(m + a + b, a, &binom2));
1432e6a796c3SToby Isaac       gra = 1. / (binom1 * binom2);
1433e6a796c3SToby Isaac     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "lgamma() - math routine is unavailable.");
1434e6a796c3SToby Isaac   }
1435e6a796c3SToby Isaac #endif
1436e6a796c3SToby Isaac   *leftw  = twoab1 * grb / b;
1437e6a796c3SToby Isaac   *rightw = twoab1 * gra / a;
14383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1439e6a796c3SToby Isaac }
1440e6a796c3SToby Isaac 
1441e6a796c3SToby Isaac /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
1442e6a796c3SToby Isaac    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
1443d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
1444d71ae5a4SJacob Faibussowitsch {
144594e21283SToby Isaac   PetscReal pn1, pn2;
144694e21283SToby Isaac   PetscReal cnm1, cnm1x, cnm2;
1447e6a796c3SToby Isaac   PetscInt  k;
1448e6a796c3SToby Isaac 
1449e6a796c3SToby Isaac   PetscFunctionBegin;
14509371c9d4SSatish Balay   if (!n) {
14519371c9d4SSatish Balay     *P = 1.0;
14523ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
14539371c9d4SSatish Balay   }
145494e21283SToby Isaac   PetscDTJacobiRecurrence_Internal(1, a, b, cnm1, cnm1x, cnm2);
145594e21283SToby Isaac   pn2 = 1.;
145694e21283SToby Isaac   pn1 = cnm1 + cnm1x * x;
14579371c9d4SSatish Balay   if (n == 1) {
14589371c9d4SSatish Balay     *P = pn1;
14593ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
14609371c9d4SSatish Balay   }
1461e6a796c3SToby Isaac   *P = 0.0;
1462e6a796c3SToby Isaac   for (k = 2; k < n + 1; ++k) {
146394e21283SToby Isaac     PetscDTJacobiRecurrence_Internal(k, a, b, cnm1, cnm1x, cnm2);
1464e6a796c3SToby Isaac 
146594e21283SToby Isaac     *P  = (cnm1 + cnm1x * x) * pn1 - cnm2 * pn2;
1466e6a796c3SToby Isaac     pn2 = pn1;
1467e6a796c3SToby Isaac     pn1 = *P;
1468e6a796c3SToby Isaac   }
14693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1470e6a796c3SToby Isaac }
1471e6a796c3SToby Isaac 
1472e6a796c3SToby Isaac /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
1473d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscInt k, PetscReal *P)
1474d71ae5a4SJacob Faibussowitsch {
1475e6a796c3SToby Isaac   PetscReal nP;
1476e6a796c3SToby Isaac   PetscInt  i;
1477e6a796c3SToby Isaac 
1478e6a796c3SToby Isaac   PetscFunctionBegin;
147917a42bb7SSatish Balay   *P = 0.0;
14803ba16761SJacob Faibussowitsch   if (k > n) PetscFunctionReturn(PETSC_SUCCESS);
14819566063dSJacob Faibussowitsch   PetscCall(PetscDTComputeJacobi(a + k, b + k, n - k, x, &nP));
1482e6a796c3SToby Isaac   for (i = 0; i < k; i++) nP *= (a + b + n + 1. + i) * 0.5;
1483e6a796c3SToby Isaac   *P = nP;
14843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1485e6a796c3SToby Isaac }
1486e6a796c3SToby Isaac 
1487d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTGaussJacobiQuadrature_Newton_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[])
1488d71ae5a4SJacob Faibussowitsch {
1489e6a796c3SToby Isaac   PetscInt  maxIter = 100;
149094e21283SToby Isaac   PetscReal eps     = PetscExpReal(0.75 * PetscLogReal(PETSC_MACHINE_EPSILON));
1491200b5abcSJed Brown   PetscReal a1, a6, gf;
1492e6a796c3SToby Isaac   PetscInt  k;
1493e6a796c3SToby Isaac 
1494e6a796c3SToby Isaac   PetscFunctionBegin;
1495e6a796c3SToby Isaac 
1496e6a796c3SToby Isaac   a1 = PetscPowReal(2.0, a + b + 1);
149794e21283SToby Isaac #if defined(PETSC_HAVE_LGAMMA)
1498200b5abcSJed Brown   {
1499200b5abcSJed Brown     PetscReal a2, a3, a4, a5;
150094e21283SToby Isaac     a2 = PetscLGamma(a + npoints + 1);
150194e21283SToby Isaac     a3 = PetscLGamma(b + npoints + 1);
150294e21283SToby Isaac     a4 = PetscLGamma(a + b + npoints + 1);
150394e21283SToby Isaac     a5 = PetscLGamma(npoints + 1);
150494e21283SToby Isaac     gf = PetscExpReal(a2 + a3 - (a4 + a5));
1505200b5abcSJed Brown   }
1506e6a796c3SToby Isaac #else
1507e6a796c3SToby Isaac   {
1508e6a796c3SToby Isaac     PetscInt ia, ib;
1509e6a796c3SToby Isaac 
1510e6a796c3SToby Isaac     ia = (PetscInt)a;
1511e6a796c3SToby Isaac     ib = (PetscInt)b;
151294e21283SToby Isaac     gf = 1.;
151394e21283SToby Isaac     if (ia == a && ia >= 0) { /* compute ratio of rising factorals wrt a */
151494e21283SToby Isaac       for (k = 0; k < ia; k++) gf *= (npoints + 1. + k) / (npoints + b + 1. + k);
151594e21283SToby Isaac     } else if (b == b && ib >= 0) { /* compute ratio of rising factorials wrt b */
151694e21283SToby Isaac       for (k = 0; k < ib; k++) gf *= (npoints + 1. + k) / (npoints + a + 1. + k);
151794e21283SToby Isaac     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "lgamma() - math routine is unavailable.");
1518e6a796c3SToby Isaac   }
1519e6a796c3SToby Isaac #endif
1520e6a796c3SToby Isaac 
152194e21283SToby Isaac   a6 = a1 * gf;
1522e6a796c3SToby Isaac   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
1523e6a796c3SToby Isaac    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
1524e6a796c3SToby Isaac   for (k = 0; k < npoints; ++k) {
152594e21283SToby Isaac     PetscReal r = PetscCosReal(PETSC_PI * (1. - (4. * k + 3. + 2. * b) / (4. * npoints + 2. * (a + b + 1.)))), dP;
1526e6a796c3SToby Isaac     PetscInt  j;
1527e6a796c3SToby Isaac 
1528e6a796c3SToby Isaac     if (k > 0) r = 0.5 * (r + x[k - 1]);
1529e6a796c3SToby Isaac     for (j = 0; j < maxIter; ++j) {
1530e6a796c3SToby Isaac       PetscReal s = 0.0, delta, f, fp;
1531e6a796c3SToby Isaac       PetscInt  i;
1532e6a796c3SToby Isaac 
1533e6a796c3SToby Isaac       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
15349566063dSJacob Faibussowitsch       PetscCall(PetscDTComputeJacobi(a, b, npoints, r, &f));
15359566063dSJacob Faibussowitsch       PetscCall(PetscDTComputeJacobiDerivative(a, b, npoints, r, 1, &fp));
1536e6a796c3SToby Isaac       delta = f / (fp - f * s);
1537e6a796c3SToby Isaac       r     = r - delta;
1538e6a796c3SToby Isaac       if (PetscAbsReal(delta) < eps) break;
1539e6a796c3SToby Isaac     }
1540e6a796c3SToby Isaac     x[k] = r;
15419566063dSJacob Faibussowitsch     PetscCall(PetscDTComputeJacobiDerivative(a, b, npoints, x[k], 1, &dP));
1542e6a796c3SToby Isaac     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
1543e6a796c3SToby Isaac   }
15443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1545e6a796c3SToby Isaac }
1546e6a796c3SToby Isaac 
154794e21283SToby Isaac /* Compute the diagonals of the Jacobi matrix used in Golub & Welsch algorithms for Gauss-Jacobi
1548e6a796c3SToby Isaac  * quadrature weight calculations on [-1,1] for exponents (1. + x)^a (1.-x)^b */
1549d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTJacobiMatrix_Internal(PetscInt nPoints, PetscReal a, PetscReal b, PetscReal *d, PetscReal *s)
1550d71ae5a4SJacob Faibussowitsch {
1551e6a796c3SToby Isaac   PetscInt i;
1552e6a796c3SToby Isaac 
1553e6a796c3SToby Isaac   PetscFunctionBegin;
1554e6a796c3SToby Isaac   for (i = 0; i < nPoints; i++) {
155594e21283SToby Isaac     PetscReal A, B, C;
1556e6a796c3SToby Isaac 
155794e21283SToby Isaac     PetscDTJacobiRecurrence_Internal(i + 1, a, b, A, B, C);
155894e21283SToby Isaac     d[i] = -A / B;
155994e21283SToby Isaac     if (i) s[i - 1] *= C / B;
156094e21283SToby Isaac     if (i < nPoints - 1) s[i] = 1. / B;
1561e6a796c3SToby Isaac   }
15623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1563e6a796c3SToby Isaac }
1564e6a796c3SToby Isaac 
1565d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
1566d71ae5a4SJacob Faibussowitsch {
1567e6a796c3SToby Isaac   PetscReal mu0;
1568e6a796c3SToby Isaac   PetscReal ga, gb, gab;
1569e6a796c3SToby Isaac   PetscInt  i;
1570e6a796c3SToby Isaac 
1571e6a796c3SToby Isaac   PetscFunctionBegin;
15729566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(GolubWelschCitation, &GolubWelschCite));
1573e6a796c3SToby Isaac 
1574e6a796c3SToby Isaac #if defined(PETSC_HAVE_TGAMMA)
1575e6a796c3SToby Isaac   ga  = PetscTGamma(a + 1);
1576e6a796c3SToby Isaac   gb  = PetscTGamma(b + 1);
1577e6a796c3SToby Isaac   gab = PetscTGamma(a + b + 2);
1578e6a796c3SToby Isaac #else
1579e6a796c3SToby Isaac   {
1580e6a796c3SToby Isaac     PetscInt ia, ib;
1581e6a796c3SToby Isaac 
1582e6a796c3SToby Isaac     ia = (PetscInt)a;
1583e6a796c3SToby Isaac     ib = (PetscInt)b;
1584e6a796c3SToby Isaac     if (ia == a && ib == b && ia + 1 > 0 && ib + 1 > 0 && ia + ib + 2 > 0) { /* All gamma(x) terms are (x-1)! terms */
15859566063dSJacob Faibussowitsch       PetscCall(PetscDTFactorial(ia, &ga));
15869566063dSJacob Faibussowitsch       PetscCall(PetscDTFactorial(ib, &gb));
15879566063dSJacob Faibussowitsch       PetscCall(PetscDTFactorial(ia + ib + 1, &gb));
1588e6a796c3SToby Isaac     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "tgamma() - math routine is unavailable.");
1589e6a796c3SToby Isaac   }
1590e6a796c3SToby Isaac #endif
1591e6a796c3SToby Isaac   mu0 = PetscPowReal(2., a + b + 1.) * ga * gb / gab;
1592e6a796c3SToby Isaac 
1593e6a796c3SToby Isaac #if defined(PETSCDTGAUSSIANQUADRATURE_EIG)
1594e6a796c3SToby Isaac   {
1595e6a796c3SToby Isaac     PetscReal   *diag, *subdiag;
1596e6a796c3SToby Isaac     PetscScalar *V;
1597e6a796c3SToby Isaac 
15989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(npoints, &diag, npoints, &subdiag));
15999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(npoints * npoints, &V));
16009566063dSJacob Faibussowitsch     PetscCall(PetscDTJacobiMatrix_Internal(npoints, a, b, diag, subdiag));
1601e6a796c3SToby Isaac     for (i = 0; i < npoints - 1; i++) subdiag[i] = PetscSqrtReal(subdiag[i]);
16029566063dSJacob Faibussowitsch     PetscCall(PetscDTSymmetricTridiagonalEigensolve(npoints, diag, subdiag, x, V));
160394e21283SToby Isaac     for (i = 0; i < npoints; i++) w[i] = PetscSqr(PetscRealPart(V[i * npoints])) * mu0;
16049566063dSJacob Faibussowitsch     PetscCall(PetscFree(V));
16059566063dSJacob Faibussowitsch     PetscCall(PetscFree2(diag, subdiag));
1606e6a796c3SToby Isaac   }
1607e6a796c3SToby Isaac #else
1608e6a796c3SToby Isaac   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
1609e6a796c3SToby Isaac #endif
161094e21283SToby Isaac   { /* As of March 2, 2020, The Sun Performance Library breaks the LAPACK contract for xstegr and xsteqr: the
161194e21283SToby Isaac        eigenvalues are not guaranteed to be in ascending order.  So we heave a passive aggressive sigh and check that
161294e21283SToby Isaac        the eigenvalues are sorted */
161394e21283SToby Isaac     PetscBool sorted;
161494e21283SToby Isaac 
16159566063dSJacob Faibussowitsch     PetscCall(PetscSortedReal(npoints, x, &sorted));
161694e21283SToby Isaac     if (!sorted) {
161794e21283SToby Isaac       PetscInt  *order, i;
161894e21283SToby Isaac       PetscReal *tmp;
161994e21283SToby Isaac 
16209566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(npoints, &order, npoints, &tmp));
162194e21283SToby Isaac       for (i = 0; i < npoints; i++) order[i] = i;
16229566063dSJacob Faibussowitsch       PetscCall(PetscSortRealWithPermutation(npoints, x, order));
16239566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(tmp, x, npoints));
162494e21283SToby Isaac       for (i = 0; i < npoints; i++) x[i] = tmp[order[i]];
16259566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(tmp, w, npoints));
162694e21283SToby Isaac       for (i = 0; i < npoints; i++) w[i] = tmp[order[i]];
16279566063dSJacob Faibussowitsch       PetscCall(PetscFree2(order, tmp));
162894e21283SToby Isaac     }
162994e21283SToby Isaac   }
16303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1631e6a796c3SToby Isaac }
1632e6a796c3SToby Isaac 
1633d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTGaussJacobiQuadrature_Internal(PetscInt npoints, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
1634d71ae5a4SJacob Faibussowitsch {
1635e6a796c3SToby Isaac   PetscFunctionBegin;
163608401ef6SPierre Jolivet   PetscCheck(npoints >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be positive");
1637e6a796c3SToby Isaac   /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
163808401ef6SPierre Jolivet   PetscCheck(alpha > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "alpha must be > -1.");
163908401ef6SPierre Jolivet   PetscCheck(beta > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "beta must be > -1.");
1640e6a796c3SToby Isaac 
16411baa6e33SBarry Smith   if (newton) PetscCall(PetscDTGaussJacobiQuadrature_Newton_Internal(npoints, alpha, beta, x, w));
16421baa6e33SBarry Smith   else PetscCall(PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(npoints, alpha, beta, x, w));
1643e6a796c3SToby Isaac   if (alpha == beta) { /* symmetrize */
1644e6a796c3SToby Isaac     PetscInt i;
1645e6a796c3SToby Isaac     for (i = 0; i < (npoints + 1) / 2; i++) {
1646e6a796c3SToby Isaac       PetscInt  j  = npoints - 1 - i;
1647e6a796c3SToby Isaac       PetscReal xi = x[i];
1648e6a796c3SToby Isaac       PetscReal xj = x[j];
1649e6a796c3SToby Isaac       PetscReal wi = w[i];
1650e6a796c3SToby Isaac       PetscReal wj = w[j];
1651e6a796c3SToby Isaac 
1652e6a796c3SToby Isaac       x[i] = (xi - xj) / 2.;
1653e6a796c3SToby Isaac       x[j] = (xj - xi) / 2.;
1654e6a796c3SToby Isaac       w[i] = w[j] = (wi + wj) / 2.;
1655e6a796c3SToby Isaac     }
1656e6a796c3SToby Isaac   }
16573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1658e6a796c3SToby Isaac }
1659e6a796c3SToby Isaac 
166094e21283SToby Isaac /*@
166194e21283SToby Isaac   PetscDTGaussJacobiQuadrature - quadrature for the interval [a, b] with the weight function
166294e21283SToby Isaac   $(x-a)^\alpha (x-b)^\beta$.
166394e21283SToby Isaac 
166420f4b53cSBarry Smith   Not Collective
166594e21283SToby Isaac 
166694e21283SToby Isaac   Input Parameters:
166794e21283SToby Isaac + npoints - the number of points in the quadrature rule
166894e21283SToby Isaac . a - the left endpoint of the interval
166994e21283SToby Isaac . b - the right endpoint of the interval
167094e21283SToby Isaac . alpha - the left exponent
167194e21283SToby Isaac - beta - the right exponent
167294e21283SToby Isaac 
167394e21283SToby Isaac   Output Parameters:
167420f4b53cSBarry Smith + x - array of length `npoints`, the locations of the quadrature points
167520f4b53cSBarry Smith - w - array of length `npoints`, the weights of the quadrature points
167694e21283SToby Isaac 
167794e21283SToby Isaac   Level: intermediate
167894e21283SToby Isaac 
1679dce8aebaSBarry Smith   Note:
1680dce8aebaSBarry Smith   This quadrature rule is exact for polynomials up to degree 2*npoints - 1.
1681dce8aebaSBarry Smith 
1682dce8aebaSBarry Smith .seealso: `PetscDTGaussQuadrature()`
168394e21283SToby Isaac @*/
1684d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt npoints, PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
1685d71ae5a4SJacob Faibussowitsch {
168694e21283SToby Isaac   PetscInt i;
1687e6a796c3SToby Isaac 
1688e6a796c3SToby Isaac   PetscFunctionBegin;
16899566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal));
169094e21283SToby Isaac   if (a != -1. || b != 1.) { /* shift */
169194e21283SToby Isaac     for (i = 0; i < npoints; i++) {
169294e21283SToby Isaac       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
169394e21283SToby Isaac       w[i] *= (b - a) / 2.;
169494e21283SToby Isaac     }
169594e21283SToby Isaac   }
16963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1697e6a796c3SToby Isaac }
1698e6a796c3SToby Isaac 
1699d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTGaussLobattoJacobiQuadrature_Internal(PetscInt npoints, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
1700d71ae5a4SJacob Faibussowitsch {
1701e6a796c3SToby Isaac   PetscInt i;
1702e6a796c3SToby Isaac 
1703e6a796c3SToby Isaac   PetscFunctionBegin;
170408401ef6SPierre Jolivet   PetscCheck(npoints >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be positive");
1705e6a796c3SToby Isaac   /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
170608401ef6SPierre Jolivet   PetscCheck(alpha > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "alpha must be > -1.");
170708401ef6SPierre Jolivet   PetscCheck(beta > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "beta must be > -1.");
1708e6a796c3SToby Isaac 
1709e6a796c3SToby Isaac   x[0]           = -1.;
1710e6a796c3SToby Isaac   x[npoints - 1] = 1.;
171148a46eb9SPierre Jolivet   if (npoints > 2) PetscCall(PetscDTGaussJacobiQuadrature_Internal(npoints - 2, alpha + 1., beta + 1., &x[1], &w[1], newton));
1712ad540459SPierre Jolivet   for (i = 1; i < npoints - 1; i++) w[i] /= (1. - x[i] * x[i]);
17139566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussLobattoJacobiEndweights_Internal(npoints, alpha, beta, &w[0], &w[npoints - 1]));
17143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1715e6a796c3SToby Isaac }
1716e6a796c3SToby Isaac 
171737045ce4SJed Brown /*@
171894e21283SToby Isaac   PetscDTGaussLobattoJacobiQuadrature - quadrature for the interval [a, b] with the weight function
171994e21283SToby Isaac   $(x-a)^\alpha (x-b)^\beta$, with endpoints a and b included as quadrature points.
172094e21283SToby Isaac 
172120f4b53cSBarry Smith   Not Collective
172294e21283SToby Isaac 
172394e21283SToby Isaac   Input Parameters:
172494e21283SToby Isaac + npoints - the number of points in the quadrature rule
172594e21283SToby Isaac . a - the left endpoint of the interval
172694e21283SToby Isaac . b - the right endpoint of the interval
172794e21283SToby Isaac . alpha - the left exponent
172894e21283SToby Isaac - beta - the right exponent
172994e21283SToby Isaac 
173094e21283SToby Isaac   Output Parameters:
173120f4b53cSBarry Smith + x - array of length `npoints`, the locations of the quadrature points
173220f4b53cSBarry Smith - w - array of length `npoints`, the weights of the quadrature points
173394e21283SToby Isaac 
173494e21283SToby Isaac   Level: intermediate
173594e21283SToby Isaac 
1736dce8aebaSBarry Smith   Note:
1737dce8aebaSBarry Smith   This quadrature rule is exact for polynomials up to degree 2*npoints - 3.
1738dce8aebaSBarry Smith 
1739dce8aebaSBarry Smith .seealso: `PetscDTGaussJacobiQuadrature()`
174094e21283SToby Isaac @*/
1741d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt npoints, PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
1742d71ae5a4SJacob Faibussowitsch {
174394e21283SToby Isaac   PetscInt i;
174494e21283SToby Isaac 
174594e21283SToby Isaac   PetscFunctionBegin;
17469566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal));
174794e21283SToby Isaac   if (a != -1. || b != 1.) { /* shift */
174894e21283SToby Isaac     for (i = 0; i < npoints; i++) {
174994e21283SToby Isaac       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
175094e21283SToby Isaac       w[i] *= (b - a) / 2.;
175194e21283SToby Isaac     }
175294e21283SToby Isaac   }
17533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
175494e21283SToby Isaac }
175594e21283SToby Isaac 
175694e21283SToby Isaac /*@
1757e6a796c3SToby Isaac    PetscDTGaussQuadrature - create Gauss-Legendre quadrature
175837045ce4SJed Brown 
175937045ce4SJed Brown    Not Collective
176037045ce4SJed Brown 
17614165533cSJose E. Roman    Input Parameters:
176237045ce4SJed Brown +  npoints - number of points
176337045ce4SJed Brown .  a - left end of interval (often-1)
176437045ce4SJed Brown -  b - right end of interval (often +1)
176537045ce4SJed Brown 
17664165533cSJose E. Roman    Output Parameters:
176737045ce4SJed Brown +  x - quadrature points
176837045ce4SJed Brown -  w - quadrature weights
176937045ce4SJed Brown 
177037045ce4SJed Brown    Level: intermediate
177137045ce4SJed Brown 
177237045ce4SJed Brown    References:
1773606c0280SSatish Balay .  * - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
177437045ce4SJed Brown 
1775dce8aebaSBarry Smith .seealso: `PetscDTLegendreEval()`, `PetscDTGaussJacobiQuadrature()`
177637045ce4SJed Brown @*/
1777d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
1778d71ae5a4SJacob Faibussowitsch {
177937045ce4SJed Brown   PetscInt i;
178037045ce4SJed Brown 
178137045ce4SJed Brown   PetscFunctionBegin;
17829566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussJacobiQuadrature_Internal(npoints, 0., 0., x, w, PetscDTGaussQuadratureNewton_Internal));
178394e21283SToby Isaac   if (a != -1. || b != 1.) { /* shift */
178437045ce4SJed Brown     for (i = 0; i < npoints; i++) {
1785e6a796c3SToby Isaac       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1786e6a796c3SToby Isaac       w[i] *= (b - a) / 2.;
178737045ce4SJed Brown     }
178837045ce4SJed Brown   }
17893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
179037045ce4SJed Brown }
1791194825f6SJed Brown 
17928272889dSSatish Balay /*@C
17938272889dSSatish Balay    PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre
17948272889dSSatish Balay                       nodes of a given size on the domain [-1,1]
17958272889dSSatish Balay 
17968272889dSSatish Balay    Not Collective
17978272889dSSatish Balay 
1798d8d19677SJose E. Roman    Input Parameters:
17998272889dSSatish Balay +  n - number of grid nodes
1800dce8aebaSBarry Smith -  type - `PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA` or `PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON`
18018272889dSSatish Balay 
18024165533cSJose E. Roman    Output Parameters:
18038272889dSSatish Balay +  x - quadrature points
18048272889dSSatish Balay -  w - quadrature weights
18058272889dSSatish Balay 
1806dce8aebaSBarry Smith    Level: intermediate
1807dce8aebaSBarry Smith 
18088272889dSSatish Balay    Notes:
18098272889dSSatish Balay     For n > 30  the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
18108272889dSSatish Balay           close enough to the desired solution
18118272889dSSatish Balay 
18128272889dSSatish Balay    These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes
18138272889dSSatish Balay 
1814a8d69d7bSBarry Smith    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
18158272889dSSatish Balay 
1816dce8aebaSBarry Smith .seealso: `PetscDTGaussQuadrature()`, `PetscGaussLobattoLegendreCreateType`
18178272889dSSatish Balay 
18188272889dSSatish Balay @*/
1819d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints, PetscGaussLobattoLegendreCreateType type, PetscReal *x, PetscReal *w)
1820d71ae5a4SJacob Faibussowitsch {
1821e6a796c3SToby Isaac   PetscBool newton;
18228272889dSSatish Balay 
18238272889dSSatish Balay   PetscFunctionBegin;
182408401ef6SPierre Jolivet   PetscCheck(npoints >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must provide at least 2 grid points per element");
182594e21283SToby Isaac   newton = (PetscBool)(type == PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON);
18269566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, 0., 0., x, w, newton));
18273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18288272889dSSatish Balay }
18298272889dSSatish Balay 
1830744bafbcSMatthew G. Knepley /*@
1831744bafbcSMatthew G. Knepley   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
1832744bafbcSMatthew G. Knepley 
1833744bafbcSMatthew G. Knepley   Not Collective
1834744bafbcSMatthew G. Knepley 
18354165533cSJose E. Roman   Input Parameters:
1836744bafbcSMatthew G. Knepley + dim     - The spatial dimension
1837a6b92713SMatthew G. Knepley . Nc      - The number of components
1838744bafbcSMatthew G. Knepley . npoints - number of points in one dimension
1839744bafbcSMatthew G. Knepley . a       - left end of interval (often-1)
1840744bafbcSMatthew G. Knepley - b       - right end of interval (often +1)
1841744bafbcSMatthew G. Knepley 
18424165533cSJose E. Roman   Output Parameter:
1843dce8aebaSBarry Smith . q - A `PetscQuadrature` object
1844744bafbcSMatthew G. Knepley 
1845744bafbcSMatthew G. Knepley   Level: intermediate
1846744bafbcSMatthew G. Knepley 
1847db781477SPatrick Sanan .seealso: `PetscDTGaussQuadrature()`, `PetscDTLegendreEval()`
1848744bafbcSMatthew G. Knepley @*/
1849d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1850d71ae5a4SJacob Faibussowitsch {
1851a6b92713SMatthew G. Knepley   PetscInt   totpoints = dim > 1 ? dim > 2 ? npoints * PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
1852744bafbcSMatthew G. Knepley   PetscReal *x, *w, *xw, *ww;
1853744bafbcSMatthew G. Knepley 
1854744bafbcSMatthew G. Knepley   PetscFunctionBegin;
18559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(totpoints * dim, &x));
18569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(totpoints * Nc, &w));
1857744bafbcSMatthew G. Knepley   /* Set up the Golub-Welsch system */
1858744bafbcSMatthew G. Knepley   switch (dim) {
1859744bafbcSMatthew G. Knepley   case 0:
18609566063dSJacob Faibussowitsch     PetscCall(PetscFree(x));
18619566063dSJacob Faibussowitsch     PetscCall(PetscFree(w));
18629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &x));
18639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nc, &w));
1864744bafbcSMatthew G. Knepley     x[0] = 0.0;
1865a6b92713SMatthew G. Knepley     for (c = 0; c < Nc; ++c) w[c] = 1.0;
1866744bafbcSMatthew G. Knepley     break;
1867744bafbcSMatthew G. Knepley   case 1:
18689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(npoints, &ww));
18699566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussQuadrature(npoints, a, b, x, ww));
18709371c9d4SSatish Balay     for (i = 0; i < npoints; ++i)
18719371c9d4SSatish Balay       for (c = 0; c < Nc; ++c) w[i * Nc + c] = ww[i];
18729566063dSJacob Faibussowitsch     PetscCall(PetscFree(ww));
1873744bafbcSMatthew G. Knepley     break;
1874744bafbcSMatthew G. Knepley   case 2:
18759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(npoints, &xw, npoints, &ww));
18769566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussQuadrature(npoints, a, b, xw, ww));
1877744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
1878744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
1879744bafbcSMatthew G. Knepley         x[(i * npoints + j) * dim + 0] = xw[i];
1880744bafbcSMatthew G. Knepley         x[(i * npoints + j) * dim + 1] = xw[j];
1881a6b92713SMatthew G. Knepley         for (c = 0; c < Nc; ++c) w[(i * npoints + j) * Nc + c] = ww[i] * ww[j];
1882744bafbcSMatthew G. Knepley       }
1883744bafbcSMatthew G. Knepley     }
18849566063dSJacob Faibussowitsch     PetscCall(PetscFree2(xw, ww));
1885744bafbcSMatthew G. Knepley     break;
1886744bafbcSMatthew G. Knepley   case 3:
18879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(npoints, &xw, npoints, &ww));
18889566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussQuadrature(npoints, a, b, xw, ww));
1889744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
1890744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
1891744bafbcSMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
1892744bafbcSMatthew G. Knepley           x[((i * npoints + j) * npoints + k) * dim + 0] = xw[i];
1893744bafbcSMatthew G. Knepley           x[((i * npoints + j) * npoints + k) * dim + 1] = xw[j];
1894744bafbcSMatthew G. Knepley           x[((i * npoints + j) * npoints + k) * dim + 2] = xw[k];
1895a6b92713SMatthew G. Knepley           for (c = 0; c < Nc; ++c) w[((i * npoints + j) * npoints + k) * Nc + c] = ww[i] * ww[j] * ww[k];
1896744bafbcSMatthew G. Knepley         }
1897744bafbcSMatthew G. Knepley       }
1898744bafbcSMatthew G. Knepley     }
18999566063dSJacob Faibussowitsch     PetscCall(PetscFree2(xw, ww));
1900744bafbcSMatthew G. Knepley     break;
1901d71ae5a4SJacob Faibussowitsch   default:
1902d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %" PetscInt_FMT, dim);
1903744bafbcSMatthew G. Knepley   }
19049566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
19059566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetOrder(*q, 2 * npoints - 1));
19069566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w));
19079566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*q, "GaussTensor"));
19083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1909744bafbcSMatthew G. Knepley }
1910744bafbcSMatthew G. Knepley 
1911f5f57ec0SBarry Smith /*@
1912e6a796c3SToby Isaac   PetscDTStroudConicalQuadrature - create Stroud conical quadrature for a simplex
1913494e7359SMatthew G. Knepley 
1914494e7359SMatthew G. Knepley   Not Collective
1915494e7359SMatthew G. Knepley 
19164165533cSJose E. Roman   Input Parameters:
1917494e7359SMatthew G. Knepley + dim     - The simplex dimension
1918a6b92713SMatthew G. Knepley . Nc      - The number of components
1919dcce0ee2SMatthew G. Knepley . npoints - The number of points in one dimension
1920494e7359SMatthew G. Knepley . a       - left end of interval (often-1)
1921494e7359SMatthew G. Knepley - b       - right end of interval (often +1)
1922494e7359SMatthew G. Knepley 
19234165533cSJose E. Roman   Output Parameter:
192420f4b53cSBarry Smith . q - A `PetscQuadrature` object
1925494e7359SMatthew G. Knepley 
1926494e7359SMatthew G. Knepley   Level: intermediate
1927494e7359SMatthew G. Knepley 
1928dce8aebaSBarry Smith   Note:
192920f4b53cSBarry Smith   For `dim` == 1, this is Gauss-Legendre quadrature
1930dce8aebaSBarry Smith 
1931494e7359SMatthew G. Knepley   References:
1932606c0280SSatish Balay . * - Karniadakis and Sherwin.  FIAT
1933494e7359SMatthew G. Knepley 
1934db781477SPatrick Sanan .seealso: `PetscDTGaussTensorQuadrature()`, `PetscDTGaussQuadrature()`
1935494e7359SMatthew G. Knepley @*/
1936d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1937d71ae5a4SJacob Faibussowitsch {
1938fbdc3dfeSToby Isaac   PetscInt   totprev, totrem;
1939fbdc3dfeSToby Isaac   PetscInt   totpoints;
1940fbdc3dfeSToby Isaac   PetscReal *p1, *w1;
1941fbdc3dfeSToby Isaac   PetscReal *x, *w;
1942fbdc3dfeSToby Isaac   PetscInt   i, j, k, l, m, pt, c;
1943494e7359SMatthew G. Knepley 
1944494e7359SMatthew G. Knepley   PetscFunctionBegin;
194508401ef6SPierre Jolivet   PetscCheck(!(a != -1.0) && !(b != 1.0), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
1946fbdc3dfeSToby Isaac   totpoints = 1;
1947fbdc3dfeSToby Isaac   for (i = 0, totpoints = 1; i < dim; i++) totpoints *= npoints;
19489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(totpoints * dim, &x));
19499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(totpoints * Nc, &w));
19509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(npoints, &p1, npoints, &w1));
1951fbdc3dfeSToby Isaac   for (i = 0; i < totpoints * Nc; i++) w[i] = 1.;
1952fbdc3dfeSToby Isaac   for (i = 0, totprev = 1, totrem = totpoints / npoints; i < dim; i++) {
1953fbdc3dfeSToby Isaac     PetscReal mul;
1954fbdc3dfeSToby Isaac 
1955fbdc3dfeSToby Isaac     mul = PetscPowReal(2., -i);
19569566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussJacobiQuadrature(npoints, -1., 1., i, 0.0, p1, w1));
1957fbdc3dfeSToby Isaac     for (pt = 0, l = 0; l < totprev; l++) {
1958fbdc3dfeSToby Isaac       for (j = 0; j < npoints; j++) {
1959fbdc3dfeSToby Isaac         for (m = 0; m < totrem; m++, pt++) {
1960fbdc3dfeSToby Isaac           for (k = 0; k < i; k++) x[pt * dim + k] = (x[pt * dim + k] + 1.) * (1. - p1[j]) * 0.5 - 1.;
1961fbdc3dfeSToby Isaac           x[pt * dim + i] = p1[j];
1962fbdc3dfeSToby Isaac           for (c = 0; c < Nc; c++) w[pt * Nc + c] *= mul * w1[j];
1963494e7359SMatthew G. Knepley         }
1964494e7359SMatthew G. Knepley       }
1965494e7359SMatthew G. Knepley     }
1966fbdc3dfeSToby Isaac     totprev *= npoints;
1967fbdc3dfeSToby Isaac     totrem /= npoints;
1968494e7359SMatthew G. Knepley   }
19699566063dSJacob Faibussowitsch   PetscCall(PetscFree2(p1, w1));
19709566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
19719566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetOrder(*q, 2 * npoints - 1));
19729566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w));
19739566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*q, "StroudConical"));
19743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1975494e7359SMatthew G. Knepley }
1976494e7359SMatthew G. Knepley 
1977d3c69ad0SToby Isaac static PetscBool MinSymTriQuadCite       = PETSC_FALSE;
19789371c9d4SSatish Balay const char       MinSymTriQuadCitation[] = "@article{WitherdenVincent2015,\n"
1979d3c69ad0SToby Isaac                                            "  title = {On the identification of symmetric quadrature rules for finite element methods},\n"
1980d3c69ad0SToby Isaac                                            "  journal = {Computers & Mathematics with Applications},\n"
1981d3c69ad0SToby Isaac                                            "  volume = {69},\n"
1982d3c69ad0SToby Isaac                                            "  number = {10},\n"
1983d3c69ad0SToby Isaac                                            "  pages = {1232-1241},\n"
1984d3c69ad0SToby Isaac                                            "  year = {2015},\n"
1985d3c69ad0SToby Isaac                                            "  issn = {0898-1221},\n"
1986d3c69ad0SToby Isaac                                            "  doi = {10.1016/j.camwa.2015.03.017},\n"
1987d3c69ad0SToby Isaac                                            "  url = {https://www.sciencedirect.com/science/article/pii/S0898122115001224},\n"
1988d3c69ad0SToby Isaac                                            "  author = {F.D. Witherden and P.E. Vincent},\n"
1989d3c69ad0SToby Isaac                                            "}\n";
1990d3c69ad0SToby Isaac 
1991d3c69ad0SToby Isaac #include "petscdttriquadrules.h"
1992d3c69ad0SToby Isaac 
1993d3c69ad0SToby Isaac static PetscBool MinSymTetQuadCite       = PETSC_FALSE;
19949371c9d4SSatish Balay const char       MinSymTetQuadCitation[] = "@article{JaskowiecSukumar2021\n"
1995d3c69ad0SToby Isaac                                            "  author = {Jaskowiec, Jan and Sukumar, N.},\n"
1996d3c69ad0SToby Isaac                                            "  title = {High-order symmetric cubature rules for tetrahedra and pyramids},\n"
1997d3c69ad0SToby Isaac                                            "  journal = {International Journal for Numerical Methods in Engineering},\n"
1998d3c69ad0SToby Isaac                                            "  volume = {122},\n"
1999d3c69ad0SToby Isaac                                            "  number = {1},\n"
2000d3c69ad0SToby Isaac                                            "  pages = {148-171},\n"
2001d3c69ad0SToby Isaac                                            "  doi = {10.1002/nme.6528},\n"
2002d3c69ad0SToby Isaac                                            "  url = {https://onlinelibrary.wiley.com/doi/abs/10.1002/nme.6528},\n"
2003d3c69ad0SToby Isaac                                            "  eprint = {https://onlinelibrary.wiley.com/doi/pdf/10.1002/nme.6528},\n"
2004d3c69ad0SToby Isaac                                            "  year = {2021}\n"
2005d3c69ad0SToby Isaac                                            "}\n";
2006d3c69ad0SToby Isaac 
2007d3c69ad0SToby Isaac #include "petscdttetquadrules.h"
2008d3c69ad0SToby Isaac 
2009d3c69ad0SToby Isaac // https://en.wikipedia.org/wiki/Partition_(number_theory)
2010d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTPartitionNumber(PetscInt n, PetscInt *p)
2011d71ae5a4SJacob Faibussowitsch {
2012d3c69ad0SToby Isaac   // sequence A000041 in the OEIS
2013d3c69ad0SToby Isaac   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};
2014d3c69ad0SToby Isaac   PetscInt       tabulated_max = PETSC_STATIC_ARRAY_LENGTH(partition) - 1;
2015d3c69ad0SToby Isaac 
2016d3c69ad0SToby Isaac   PetscFunctionBegin;
2017d3c69ad0SToby Isaac   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Partition number not defined for negative number %" PetscInt_FMT, n);
2018d3c69ad0SToby Isaac   // not implementing the pentagonal number recurrence, we don't need partition numbers for n that high
2019d3c69ad0SToby Isaac   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);
2020d3c69ad0SToby Isaac   *p = partition[n];
20213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2022d3c69ad0SToby Isaac }
2023d3c69ad0SToby Isaac 
2024d3c69ad0SToby Isaac /*@
2025d3c69ad0SToby Isaac   PetscDTSimplexQuadrature - Create a quadrature rule for a simplex that exactly integrates polynomials up to a given degree.
2026d3c69ad0SToby Isaac 
2027d3c69ad0SToby Isaac   Not Collective
2028d3c69ad0SToby Isaac 
2029d3c69ad0SToby Isaac   Input Parameters:
2030d3c69ad0SToby Isaac + dim     - The spatial dimension of the simplex (1 = segment, 2 = triangle, 3 = tetrahedron)
2031d3c69ad0SToby Isaac . degree  - The largest polynomial degree that is required to be integrated exactly
2032d3c69ad0SToby Isaac - type    - left end of interval (often-1)
2033d3c69ad0SToby Isaac 
2034d3c69ad0SToby Isaac   Output Parameter:
2035dce8aebaSBarry Smith . quad    - A `PetscQuadrature` object for integration over the biunit simplex
2036d3c69ad0SToby Isaac             (defined by the bounds $x_i >= -1$ and $\sum_i x_i <= 2 - d$) that is exact for
2037d3c69ad0SToby Isaac             polynomials up to the given degree
2038d3c69ad0SToby Isaac 
2039d3c69ad0SToby Isaac   Level: intermediate
2040d3c69ad0SToby Isaac 
2041dce8aebaSBarry Smith .seealso: `PetscDTSimplexQuadratureType`, `PetscDTGaussQuadrature()`, `PetscDTStroudCononicalQuadrature()`, `PetscQuadrature`
2042d3c69ad0SToby Isaac @*/
2043d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTSimplexQuadrature(PetscInt dim, PetscInt degree, PetscDTSimplexQuadratureType type, PetscQuadrature *quad)
2044d71ae5a4SJacob Faibussowitsch {
2045d3c69ad0SToby Isaac   PetscDTSimplexQuadratureType orig_type = type;
2046d3c69ad0SToby Isaac 
2047d3c69ad0SToby Isaac   PetscFunctionBegin;
2048d3c69ad0SToby Isaac   PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative dimension %" PetscInt_FMT, dim);
2049d3c69ad0SToby Isaac   PetscCheck(degree >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %" PetscInt_FMT, degree);
2050ad540459SPierre Jolivet   if (type == PETSCDTSIMPLEXQUAD_DEFAULT) type = PETSCDTSIMPLEXQUAD_MINSYM;
2051d3c69ad0SToby Isaac   if (type == PETSCDTSIMPLEXQUAD_CONIC || dim < 2) {
2052d3c69ad0SToby Isaac     PetscInt points_per_dim = (degree + 2) / 2; // ceil((degree + 1) / 2);
2053d3c69ad0SToby Isaac     PetscCall(PetscDTStroudConicalQuadrature(dim, 1, points_per_dim, -1, 1, quad));
2054d3c69ad0SToby Isaac   } else {
2055d3c69ad0SToby Isaac     PetscInt          n    = dim + 1;
2056d3c69ad0SToby Isaac     PetscInt          fact = 1;
2057d3c69ad0SToby Isaac     PetscInt         *part, *perm;
2058d3c69ad0SToby Isaac     PetscInt          p = 0;
2059d3c69ad0SToby Isaac     PetscInt          max_degree;
2060d3c69ad0SToby Isaac     const PetscInt   *nodes_per_type     = NULL;
2061d3c69ad0SToby Isaac     const PetscInt   *all_num_full_nodes = NULL;
2062d3c69ad0SToby Isaac     const PetscReal **weights_list       = NULL;
2063d3c69ad0SToby Isaac     const PetscReal **compact_nodes_list = NULL;
2064d3c69ad0SToby Isaac     const char       *citation           = NULL;
2065d3c69ad0SToby Isaac     PetscBool        *cited              = NULL;
2066d3c69ad0SToby Isaac 
2067d3c69ad0SToby Isaac     switch (dim) {
2068d3c69ad0SToby Isaac     case 2:
2069d3c69ad0SToby Isaac       cited              = &MinSymTriQuadCite;
2070d3c69ad0SToby Isaac       citation           = MinSymTriQuadCitation;
2071d3c69ad0SToby Isaac       max_degree         = PetscDTWVTriQuad_max_degree;
2072d3c69ad0SToby Isaac       nodes_per_type     = PetscDTWVTriQuad_num_orbits;
2073d3c69ad0SToby Isaac       all_num_full_nodes = PetscDTWVTriQuad_num_nodes;
2074d3c69ad0SToby Isaac       weights_list       = PetscDTWVTriQuad_weights;
2075d3c69ad0SToby Isaac       compact_nodes_list = PetscDTWVTriQuad_orbits;
2076d3c69ad0SToby Isaac       break;
2077d3c69ad0SToby Isaac     case 3:
2078d3c69ad0SToby Isaac       cited              = &MinSymTetQuadCite;
2079d3c69ad0SToby Isaac       citation           = MinSymTetQuadCitation;
2080d3c69ad0SToby Isaac       max_degree         = PetscDTJSTetQuad_max_degree;
2081d3c69ad0SToby Isaac       nodes_per_type     = PetscDTJSTetQuad_num_orbits;
2082d3c69ad0SToby Isaac       all_num_full_nodes = PetscDTJSTetQuad_num_nodes;
2083d3c69ad0SToby Isaac       weights_list       = PetscDTJSTetQuad_weights;
2084d3c69ad0SToby Isaac       compact_nodes_list = PetscDTJSTetQuad_orbits;
2085d3c69ad0SToby Isaac       break;
2086d71ae5a4SJacob Faibussowitsch     default:
2087d71ae5a4SJacob Faibussowitsch       max_degree = -1;
2088d71ae5a4SJacob Faibussowitsch       break;
2089d3c69ad0SToby Isaac     }
2090d3c69ad0SToby Isaac 
2091d3c69ad0SToby Isaac     if (degree > max_degree) {
2092d3c69ad0SToby Isaac       if (orig_type == PETSCDTSIMPLEXQUAD_DEFAULT) {
2093d3c69ad0SToby Isaac         // fall back to conic
2094d3c69ad0SToby Isaac         PetscCall(PetscDTSimplexQuadrature(dim, degree, PETSCDTSIMPLEXQUAD_CONIC, quad));
20953ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
2096d3c69ad0SToby Isaac       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Minimal symmetric quadrature for dim %" PetscInt_FMT ", degree %" PetscInt_FMT " unsupported", dim, degree);
2097d3c69ad0SToby Isaac     }
2098d3c69ad0SToby Isaac 
2099d3c69ad0SToby Isaac     PetscCall(PetscCitationsRegister(citation, cited));
2100d3c69ad0SToby Isaac 
2101d3c69ad0SToby Isaac     PetscCall(PetscDTPartitionNumber(n, &p));
2102d3c69ad0SToby Isaac     for (PetscInt d = 2; d <= n; d++) fact *= d;
2103d3c69ad0SToby Isaac 
2104d3c69ad0SToby Isaac     PetscInt         num_full_nodes      = all_num_full_nodes[degree];
2105d3c69ad0SToby Isaac     const PetscReal *all_compact_nodes   = compact_nodes_list[degree];
2106d3c69ad0SToby Isaac     const PetscReal *all_compact_weights = weights_list[degree];
2107d3c69ad0SToby Isaac     nodes_per_type                       = &nodes_per_type[p * degree];
2108d3c69ad0SToby Isaac 
2109d3c69ad0SToby Isaac     PetscReal      *points;
2110d3c69ad0SToby Isaac     PetscReal      *counts;
2111d3c69ad0SToby Isaac     PetscReal      *weights;
2112d3c69ad0SToby Isaac     PetscReal      *bary_to_biunit; // row-major transformation of barycentric coordinate to biunit
2113d3c69ad0SToby Isaac     PetscQuadrature q;
2114d3c69ad0SToby Isaac 
2115d3c69ad0SToby Isaac     // compute the transformation
2116d3c69ad0SToby Isaac     PetscCall(PetscMalloc1(n * dim, &bary_to_biunit));
2117d3c69ad0SToby Isaac     for (PetscInt d = 0; d < dim; d++) {
2118ad540459SPierre Jolivet       for (PetscInt b = 0; b < n; b++) bary_to_biunit[d * n + b] = (d == b) ? 1.0 : -1.0;
2119d3c69ad0SToby Isaac     }
2120d3c69ad0SToby Isaac 
2121d3c69ad0SToby Isaac     PetscCall(PetscMalloc3(n, &part, n, &perm, n, &counts));
2122d3c69ad0SToby Isaac     PetscCall(PetscCalloc1(num_full_nodes * dim, &points));
2123d3c69ad0SToby Isaac     PetscCall(PetscMalloc1(num_full_nodes, &weights));
2124d3c69ad0SToby Isaac 
2125d3c69ad0SToby Isaac     // (0, 0, ...) is the first partition lexicographically
2126d3c69ad0SToby Isaac     PetscCall(PetscArrayzero(part, n));
2127d3c69ad0SToby Isaac     PetscCall(PetscArrayzero(counts, n));
2128d3c69ad0SToby Isaac     counts[0] = n;
2129d3c69ad0SToby Isaac 
2130d3c69ad0SToby Isaac     // for each partition
2131d3c69ad0SToby Isaac     for (PetscInt s = 0, node_offset = 0; s < p; s++) {
2132d3c69ad0SToby Isaac       PetscInt num_compact_coords = part[n - 1] + 1;
2133d3c69ad0SToby Isaac 
2134d3c69ad0SToby Isaac       const PetscReal *compact_nodes   = all_compact_nodes;
2135d3c69ad0SToby Isaac       const PetscReal *compact_weights = all_compact_weights;
2136d3c69ad0SToby Isaac       all_compact_nodes += num_compact_coords * nodes_per_type[s];
2137d3c69ad0SToby Isaac       all_compact_weights += nodes_per_type[s];
2138d3c69ad0SToby Isaac 
2139d3c69ad0SToby Isaac       // for every permutation of the vertices
2140d3c69ad0SToby Isaac       for (PetscInt f = 0; f < fact; f++) {
2141d3c69ad0SToby Isaac         PetscCall(PetscDTEnumPerm(n, f, perm, NULL));
2142d3c69ad0SToby Isaac 
2143d3c69ad0SToby Isaac         // check if it is a valid permutation
2144d3c69ad0SToby Isaac         PetscInt digit;
2145d3c69ad0SToby Isaac         for (digit = 1; digit < n; digit++) {
2146d3c69ad0SToby Isaac           // skip permutations that would duplicate a node because it has a smaller symmetry group
2147d3c69ad0SToby Isaac           if (part[digit - 1] == part[digit] && perm[digit - 1] > perm[digit]) break;
2148d3c69ad0SToby Isaac         }
2149d3c69ad0SToby Isaac         if (digit < n) continue;
2150d3c69ad0SToby Isaac 
2151d3c69ad0SToby Isaac         // create full nodes from this permutation of the compact nodes
2152d3c69ad0SToby Isaac         PetscReal *full_nodes   = &points[node_offset * dim];
2153d3c69ad0SToby Isaac         PetscReal *full_weights = &weights[node_offset];
2154d3c69ad0SToby Isaac 
2155d3c69ad0SToby Isaac         PetscCall(PetscArraycpy(full_weights, compact_weights, nodes_per_type[s]));
2156d3c69ad0SToby Isaac         for (PetscInt b = 0; b < n; b++) {
2157d3c69ad0SToby Isaac           for (PetscInt d = 0; d < dim; d++) {
2158ad540459SPierre Jolivet             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]];
2159d3c69ad0SToby Isaac           }
2160d3c69ad0SToby Isaac         }
2161d3c69ad0SToby Isaac         node_offset += nodes_per_type[s];
2162d3c69ad0SToby Isaac       }
2163d3c69ad0SToby Isaac 
2164d3c69ad0SToby Isaac       if (s < p - 1) { // Generate the next partition
2165d3c69ad0SToby Isaac         /* A partition is described by the number of coordinates that are in
2166d3c69ad0SToby Isaac          * each set of duplicates (counts) and redundantly by mapping each
2167d3c69ad0SToby Isaac          * index to its set of duplicates (part)
2168d3c69ad0SToby Isaac          *
2169d3c69ad0SToby Isaac          * Counts should always be in nonincreasing order
2170d3c69ad0SToby Isaac          *
2171d3c69ad0SToby Isaac          * We want to generate the partitions lexically by part, which means
2172d3c69ad0SToby Isaac          * finding the last index where count > 1 and reducing by 1.
2173d3c69ad0SToby Isaac          *
2174d3c69ad0SToby Isaac          * For the new counts beyond that index, we eagerly assign the remaining
2175d3c69ad0SToby Isaac          * capacity of the partition to smaller indices (ensures lexical ordering),
2176d3c69ad0SToby Isaac          * while respecting the nonincreasing invariant of the counts
2177d3c69ad0SToby Isaac          */
2178d3c69ad0SToby Isaac         PetscInt last_digit            = part[n - 1];
2179d3c69ad0SToby Isaac         PetscInt last_digit_with_extra = last_digit;
2180d3c69ad0SToby Isaac         while (counts[last_digit_with_extra] == 1) last_digit_with_extra--;
2181d3c69ad0SToby Isaac         PetscInt limit               = --counts[last_digit_with_extra];
2182d3c69ad0SToby Isaac         PetscInt total_to_distribute = last_digit - last_digit_with_extra + 1;
2183d3c69ad0SToby Isaac         for (PetscInt digit = last_digit_with_extra + 1; digit < n; digit++) {
2184d3c69ad0SToby Isaac           counts[digit] = PetscMin(limit, total_to_distribute);
2185d3c69ad0SToby Isaac           total_to_distribute -= PetscMin(limit, total_to_distribute);
2186d3c69ad0SToby Isaac         }
2187d3c69ad0SToby Isaac         for (PetscInt digit = 0, offset = 0; digit < n; digit++) {
2188d3c69ad0SToby Isaac           PetscInt count = counts[digit];
2189ad540459SPierre Jolivet           for (PetscInt c = 0; c < count; c++) part[offset++] = digit;
2190d3c69ad0SToby Isaac         }
2191d3c69ad0SToby Isaac       }
2192d3c69ad0SToby Isaac     }
2193d3c69ad0SToby Isaac     PetscCall(PetscFree3(part, perm, counts));
2194d3c69ad0SToby Isaac     PetscCall(PetscFree(bary_to_biunit));
2195d3c69ad0SToby Isaac     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &q));
2196b414c505SJed Brown     PetscCall(PetscQuadratureSetOrder(q, degree));
2197d3c69ad0SToby Isaac     PetscCall(PetscQuadratureSetData(q, dim, 1, num_full_nodes, points, weights));
2198d3c69ad0SToby Isaac     *quad = q;
2199d3c69ad0SToby Isaac   }
22003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2201d3c69ad0SToby Isaac }
2202d3c69ad0SToby Isaac 
2203f5f57ec0SBarry Smith /*@
2204b3c0f97bSTom Klotz   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
2205b3c0f97bSTom Klotz 
2206b3c0f97bSTom Klotz   Not Collective
2207b3c0f97bSTom Klotz 
22084165533cSJose E. Roman   Input Parameters:
2209b3c0f97bSTom Klotz + dim   - The cell dimension
2210b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l
2211b3c0f97bSTom Klotz . a     - left end of interval (often-1)
2212b3c0f97bSTom Klotz - b     - right end of interval (often +1)
2213b3c0f97bSTom Klotz 
22144165533cSJose E. Roman   Output Parameter:
2215dce8aebaSBarry Smith . q - A `PetscQuadrature` object
2216b3c0f97bSTom Klotz 
2217b3c0f97bSTom Klotz   Level: intermediate
2218b3c0f97bSTom Klotz 
2219dce8aebaSBarry Smith .seealso: `PetscDTGaussTensorQuadrature()`, `PetscQuadrature`
2220b3c0f97bSTom Klotz @*/
2221d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
2222d71ae5a4SJacob Faibussowitsch {
2223b3c0f97bSTom Klotz   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
2224b3c0f97bSTom Klotz   const PetscReal alpha = (b - a) / 2.;              /* Half-width of the integration interval */
2225b3c0f97bSTom Klotz   const PetscReal beta  = (b + a) / 2.;              /* Center of the integration interval */
2226b3c0f97bSTom Klotz   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
2227d84b4d08SMatthew G. Knepley   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
2228b3c0f97bSTom Klotz   PetscReal       wk = 0.5 * PETSC_PI;               /* Quadrature weight at x_k */
2229b3c0f97bSTom Klotz   PetscReal      *x, *w;
2230b3c0f97bSTom Klotz   PetscInt        K, k, npoints;
2231b3c0f97bSTom Klotz 
2232b3c0f97bSTom Klotz   PetscFunctionBegin;
223363a3b9bcSJacob Faibussowitsch   PetscCheck(dim <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not yet implemented", dim);
223428b400f6SJacob Faibussowitsch   PetscCheck(level, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
2235b3c0f97bSTom Klotz   /* Find K such that the weights are < 32 digits of precision */
2236ad540459SPierre Jolivet   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)));
22379566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
22389566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetOrder(*q, 2 * K + 1));
2239b3c0f97bSTom Klotz   npoints = 2 * K - 1;
22409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints * dim, &x));
22419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &w));
2242b3c0f97bSTom Klotz   /* Center term */
2243b3c0f97bSTom Klotz   x[0] = beta;
2244b3c0f97bSTom Klotz   w[0] = 0.5 * alpha * PETSC_PI;
2245b3c0f97bSTom Klotz   for (k = 1; k < K; ++k) {
22469add2064SThomas Klotz     wk           = 0.5 * alpha * h * PETSC_PI * PetscCoshReal(k * h) / PetscSqr(PetscCoshReal(0.5 * PETSC_PI * PetscSinhReal(k * h)));
22471118d4bcSLisandro Dalcin     xk           = PetscTanhReal(0.5 * PETSC_PI * PetscSinhReal(k * h));
2248b3c0f97bSTom Klotz     x[2 * k - 1] = -alpha * xk + beta;
2249b3c0f97bSTom Klotz     w[2 * k - 1] = wk;
2250b3c0f97bSTom Klotz     x[2 * k + 0] = alpha * xk + beta;
2251b3c0f97bSTom Klotz     w[2 * k + 0] = wk;
2252b3c0f97bSTom Klotz   }
22539566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*q, dim, 1, npoints, x, w));
22543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2255b3c0f97bSTom Klotz }
2256b3c0f97bSTom Klotz 
2257d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(const PetscReal[], void *, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, void *ctx, PetscReal *sol)
2258d71ae5a4SJacob Faibussowitsch {
2259b3c0f97bSTom Klotz   const PetscInt  p     = 16;           /* Digits of precision in the evaluation */
2260b3c0f97bSTom Klotz   const PetscReal alpha = (b - a) / 2.; /* Half-width of the integration interval */
2261b3c0f97bSTom Klotz   const PetscReal beta  = (b + a) / 2.; /* Center of the integration interval */
2262b3c0f97bSTom Klotz   PetscReal       h     = 1.0;          /* Step size, length between x_k */
2263b3c0f97bSTom Klotz   PetscInt        l     = 0;            /* Level of refinement, h = 2^{-l} */
2264b3c0f97bSTom Klotz   PetscReal       osum  = 0.0;          /* Integral on last level */
2265b3c0f97bSTom Klotz   PetscReal       psum  = 0.0;          /* Integral on the level before the last level */
2266b3c0f97bSTom Klotz   PetscReal       sum;                  /* Integral on current level */
2267446c295cSMatthew G. Knepley   PetscReal       yk;                   /* Quadrature point 1 - x_k on reference domain [-1, 1] */
2268b3c0f97bSTom Klotz   PetscReal       lx, rx;               /* Quadrature points to the left and right of 0 on the real domain [a, b] */
2269b3c0f97bSTom Klotz   PetscReal       wk;                   /* Quadrature weight at x_k */
2270b3c0f97bSTom Klotz   PetscReal       lval, rval;           /* Terms in the quadature sum to the left and right of 0 */
2271b3c0f97bSTom Klotz   PetscInt        d;                    /* Digits of precision in the integral */
2272b3c0f97bSTom Klotz 
2273b3c0f97bSTom Klotz   PetscFunctionBegin;
227408401ef6SPierre Jolivet   PetscCheck(digits > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
2275b3c0f97bSTom Klotz   /* Center term */
2276d6685f55SMatthew G. Knepley   func(&beta, ctx, &lval);
2277b3c0f97bSTom Klotz   sum = 0.5 * alpha * PETSC_PI * lval;
2278b3c0f97bSTom Klotz   /* */
2279b3c0f97bSTom Klotz   do {
2280b3c0f97bSTom Klotz     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
2281b3c0f97bSTom Klotz     PetscInt  k = 1;
2282b3c0f97bSTom Klotz 
2283b3c0f97bSTom Klotz     ++l;
228463a3b9bcSJacob Faibussowitsch     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %" PetscInt_FMT " sum: %15.15f\n", l, sum); */
2285b3c0f97bSTom Klotz     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
2286b3c0f97bSTom Klotz     psum = osum;
2287b3c0f97bSTom Klotz     osum = sum;
2288b3c0f97bSTom Klotz     h *= 0.5;
2289b3c0f97bSTom Klotz     sum *= 0.5;
2290b3c0f97bSTom Klotz     do {
22919add2064SThomas Klotz       wk = 0.5 * h * PETSC_PI * PetscCoshReal(k * h) / PetscSqr(PetscCoshReal(0.5 * PETSC_PI * PetscSinhReal(k * h)));
2292446c295cSMatthew G. Knepley       yk = 1.0 / (PetscExpReal(0.5 * PETSC_PI * PetscSinhReal(k * h)) * PetscCoshReal(0.5 * PETSC_PI * PetscSinhReal(k * h)));
2293446c295cSMatthew G. Knepley       lx = -alpha * (1.0 - yk) + beta;
2294446c295cSMatthew G. Knepley       rx = alpha * (1.0 - yk) + beta;
2295d6685f55SMatthew G. Knepley       func(&lx, ctx, &lval);
2296d6685f55SMatthew G. Knepley       func(&rx, ctx, &rval);
2297b3c0f97bSTom Klotz       lterm   = alpha * wk * lval;
2298b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
2299b3c0f97bSTom Klotz       sum += lterm;
2300b3c0f97bSTom Klotz       rterm   = alpha * wk * rval;
2301b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
2302b3c0f97bSTom Klotz       sum += rterm;
2303b3c0f97bSTom Klotz       ++k;
2304b3c0f97bSTom Klotz       /* Only need to evaluate every other point on refined levels */
2305b3c0f97bSTom Klotz       if (l != 1) ++k;
23069add2064SThomas Klotz     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
2307b3c0f97bSTom Klotz 
2308b3c0f97bSTom Klotz     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
2309b3c0f97bSTom Klotz     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
2310b3c0f97bSTom Klotz     d3 = PetscLog10Real(maxTerm) - p;
231109d48545SBarry Smith     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
231209d48545SBarry Smith     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
2313b3c0f97bSTom Klotz     d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1) / d2, 2 * d1), d3), d4)));
23149add2064SThomas Klotz   } while (d < digits && l < 12);
2315b3c0f97bSTom Klotz   *sol = sum;
2316e510cb1fSThomas Klotz 
23173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2318b3c0f97bSTom Klotz }
2319b3c0f97bSTom Klotz 
2320497880caSRichard Tran Mills #if defined(PETSC_HAVE_MPFR)
2321d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(const PetscReal[], void *, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, void *ctx, PetscReal *sol)
2322d71ae5a4SJacob Faibussowitsch {
2323e510cb1fSThomas Klotz   const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */
232429f144ccSMatthew G. Knepley   PetscInt       l            = 0; /* Level of refinement, h = 2^{-l} */
232529f144ccSMatthew G. Knepley   mpfr_t         alpha;            /* Half-width of the integration interval */
232629f144ccSMatthew G. Knepley   mpfr_t         beta;             /* Center of the integration interval */
232729f144ccSMatthew G. Knepley   mpfr_t         h;                /* Step size, length between x_k */
232829f144ccSMatthew G. Knepley   mpfr_t         osum;             /* Integral on last level */
232929f144ccSMatthew G. Knepley   mpfr_t         psum;             /* Integral on the level before the last level */
233029f144ccSMatthew G. Knepley   mpfr_t         sum;              /* Integral on current level */
233129f144ccSMatthew G. Knepley   mpfr_t         yk;               /* Quadrature point 1 - x_k on reference domain [-1, 1] */
233229f144ccSMatthew G. Knepley   mpfr_t         lx, rx;           /* Quadrature points to the left and right of 0 on the real domain [a, b] */
233329f144ccSMatthew G. Knepley   mpfr_t         wk;               /* Quadrature weight at x_k */
23341fbc92bbSMatthew G. Knepley   PetscReal      lval, rval, rtmp; /* Terms in the quadature sum to the left and right of 0 */
233529f144ccSMatthew G. Knepley   PetscInt       d;                /* Digits of precision in the integral */
233629f144ccSMatthew G. Knepley   mpfr_t         pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
233729f144ccSMatthew G. Knepley 
233829f144ccSMatthew G. Knepley   PetscFunctionBegin;
233908401ef6SPierre Jolivet   PetscCheck(digits > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
234029f144ccSMatthew G. Knepley   /* Create high precision storage */
2341c9f744b5SMatthew G. Knepley   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);
234229f144ccSMatthew G. Knepley   /* Initialization */
234329f144ccSMatthew G. Knepley   mpfr_set_d(alpha, 0.5 * (b - a), MPFR_RNDN);
234429f144ccSMatthew G. Knepley   mpfr_set_d(beta, 0.5 * (b + a), MPFR_RNDN);
234529f144ccSMatthew G. Knepley   mpfr_set_d(osum, 0.0, MPFR_RNDN);
234629f144ccSMatthew G. Knepley   mpfr_set_d(psum, 0.0, MPFR_RNDN);
234729f144ccSMatthew G. Knepley   mpfr_set_d(h, 1.0, MPFR_RNDN);
234829f144ccSMatthew G. Knepley   mpfr_const_pi(pi2, MPFR_RNDN);
234929f144ccSMatthew G. Knepley   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
235029f144ccSMatthew G. Knepley   /* Center term */
23511fbc92bbSMatthew G. Knepley   rtmp = 0.5 * (b + a);
23521fbc92bbSMatthew G. Knepley   func(&rtmp, ctx, &lval);
235329f144ccSMatthew G. Knepley   mpfr_set(sum, pi2, MPFR_RNDN);
235429f144ccSMatthew G. Knepley   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
235529f144ccSMatthew G. Knepley   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
235629f144ccSMatthew G. Knepley   /* */
235729f144ccSMatthew G. Knepley   do {
235829f144ccSMatthew G. Knepley     PetscReal d1, d2, d3, d4;
235929f144ccSMatthew G. Knepley     PetscInt  k = 1;
236029f144ccSMatthew G. Knepley 
236129f144ccSMatthew G. Knepley     ++l;
236229f144ccSMatthew G. Knepley     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
236363a3b9bcSJacob Faibussowitsch     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %" PetscInt_FMT " sum: %15.15f\n", l, sum); */
236429f144ccSMatthew G. Knepley     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
236529f144ccSMatthew G. Knepley     mpfr_set(psum, osum, MPFR_RNDN);
236629f144ccSMatthew G. Knepley     mpfr_set(osum, sum, MPFR_RNDN);
236729f144ccSMatthew G. Knepley     mpfr_mul_d(h, h, 0.5, MPFR_RNDN);
236829f144ccSMatthew G. Knepley     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
236929f144ccSMatthew G. Knepley     do {
237029f144ccSMatthew G. Knepley       mpfr_set_si(kh, k, MPFR_RNDN);
237129f144ccSMatthew G. Knepley       mpfr_mul(kh, kh, h, MPFR_RNDN);
237229f144ccSMatthew G. Knepley       /* Weight */
237329f144ccSMatthew G. Knepley       mpfr_set(wk, h, MPFR_RNDN);
237429f144ccSMatthew G. Knepley       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
237529f144ccSMatthew G. Knepley       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
237629f144ccSMatthew G. Knepley       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
237729f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
237829f144ccSMatthew G. Knepley       mpfr_sqr(tmp, tmp, MPFR_RNDN);
237929f144ccSMatthew G. Knepley       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
238029f144ccSMatthew G. Knepley       mpfr_div(wk, wk, tmp, MPFR_RNDN);
238129f144ccSMatthew G. Knepley       /* Abscissa */
238229f144ccSMatthew G. Knepley       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
238329f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
238429f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
238529f144ccSMatthew G. Knepley       mpfr_exp(tmp, msinh, MPFR_RNDN);
238629f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
238729f144ccSMatthew G. Knepley       /* Quadrature points */
238829f144ccSMatthew G. Knepley       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
238929f144ccSMatthew G. Knepley       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
239029f144ccSMatthew G. Knepley       mpfr_add(lx, lx, beta, MPFR_RNDU);
239129f144ccSMatthew G. Knepley       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
239229f144ccSMatthew G. Knepley       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
239329f144ccSMatthew G. Knepley       mpfr_add(rx, rx, beta, MPFR_RNDD);
239429f144ccSMatthew G. Knepley       /* Evaluation */
23951fbc92bbSMatthew G. Knepley       rtmp = mpfr_get_d(lx, MPFR_RNDU);
23961fbc92bbSMatthew G. Knepley       func(&rtmp, ctx, &lval);
23971fbc92bbSMatthew G. Knepley       rtmp = mpfr_get_d(rx, MPFR_RNDD);
23981fbc92bbSMatthew G. Knepley       func(&rtmp, ctx, &rval);
239929f144ccSMatthew G. Knepley       /* Update */
240029f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
240129f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
240229f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
240329f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
240429f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
240529f144ccSMatthew G. Knepley       mpfr_set(curTerm, tmp, MPFR_RNDN);
240629f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
240729f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
240829f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
240929f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
241029f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
241129f144ccSMatthew G. Knepley       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
241229f144ccSMatthew G. Knepley       ++k;
241329f144ccSMatthew G. Knepley       /* Only need to evaluate every other point on refined levels */
241429f144ccSMatthew G. Knepley       if (l != 1) ++k;
241529f144ccSMatthew G. Knepley       mpfr_log10(tmp, wk, MPFR_RNDN);
241629f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
2417c9f744b5SMatthew G. Knepley     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor * digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
241829f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
241929f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
242029f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
242129f144ccSMatthew G. Knepley     d1 = mpfr_get_d(tmp, MPFR_RNDN);
242229f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
242329f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
242429f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
242529f144ccSMatthew G. Knepley     d2 = mpfr_get_d(tmp, MPFR_RNDN);
242629f144ccSMatthew G. Knepley     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
2427c9f744b5SMatthew G. Knepley     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
242829f144ccSMatthew G. Knepley     mpfr_log10(tmp, curTerm, MPFR_RNDN);
242929f144ccSMatthew G. Knepley     d4 = mpfr_get_d(tmp, MPFR_RNDN);
243029f144ccSMatthew G. Knepley     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1) / d2, 2 * d1), d3), d4)));
2431b0649871SThomas Klotz   } while (d < digits && l < 8);
243229f144ccSMatthew G. Knepley   *sol = mpfr_get_d(sum, MPFR_RNDN);
243329f144ccSMatthew G. Knepley   /* Cleanup */
243429f144ccSMatthew G. Knepley   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
24353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243629f144ccSMatthew G. Knepley }
2437d525116cSMatthew G. Knepley #else
2438fbfcfee5SBarry Smith 
2439d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(const PetscReal[], void *, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, void *ctx, PetscReal *sol)
2440d71ae5a4SJacob Faibussowitsch {
2441d525116cSMatthew G. Knepley   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
2442d525116cSMatthew G. Knepley }
244329f144ccSMatthew G. Knepley #endif
244429f144ccSMatthew G. Knepley 
24452df84da0SMatthew G. Knepley /*@
24462df84da0SMatthew G. Knepley   PetscDTTensorQuadratureCreate - create the tensor product quadrature from two lower-dimensional quadratures
24472df84da0SMatthew G. Knepley 
24482df84da0SMatthew G. Knepley   Not Collective
24492df84da0SMatthew G. Knepley 
24502df84da0SMatthew G. Knepley   Input Parameters:
24512df84da0SMatthew G. Knepley + q1 - The first quadrature
24522df84da0SMatthew G. Knepley - q2 - The second quadrature
24532df84da0SMatthew G. Knepley 
24542df84da0SMatthew G. Knepley   Output Parameter:
2455dce8aebaSBarry Smith . q - A `PetscQuadrature` object
24562df84da0SMatthew G. Knepley 
24572df84da0SMatthew G. Knepley   Level: intermediate
24582df84da0SMatthew G. Knepley 
2459dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscDTGaussTensorQuadrature()`
24602df84da0SMatthew G. Knepley @*/
2461d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTTensorQuadratureCreate(PetscQuadrature q1, PetscQuadrature q2, PetscQuadrature *q)
2462d71ae5a4SJacob Faibussowitsch {
24632df84da0SMatthew G. Knepley   const PetscReal *x1, *w1, *x2, *w2;
24642df84da0SMatthew G. Knepley   PetscReal       *x, *w;
24652df84da0SMatthew G. Knepley   PetscInt         dim1, Nc1, Np1, order1, qa, d1;
24662df84da0SMatthew G. Knepley   PetscInt         dim2, Nc2, Np2, order2, qb, d2;
24672df84da0SMatthew G. Knepley   PetscInt         dim, Nc, Np, order, qc, d;
24682df84da0SMatthew G. Knepley 
24692df84da0SMatthew G. Knepley   PetscFunctionBegin;
24702df84da0SMatthew G. Knepley   PetscValidHeaderSpecific(q1, PETSCQUADRATURE_CLASSID, 1);
24712df84da0SMatthew G. Knepley   PetscValidHeaderSpecific(q2, PETSCQUADRATURE_CLASSID, 2);
24722df84da0SMatthew G. Knepley   PetscValidPointer(q, 3);
24739566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetOrder(q1, &order1));
24749566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetOrder(q2, &order2));
24752df84da0SMatthew G. Knepley   PetscCheck(order1 == order2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Order1 %" PetscInt_FMT " != %" PetscInt_FMT " Order2", order1, order2);
24769566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(q1, &dim1, &Nc1, &Np1, &x1, &w1));
24779566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(q2, &dim2, &Nc2, &Np2, &x2, &w2));
24782df84da0SMatthew G. Knepley   PetscCheck(Nc1 == Nc2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "NumComp1 %" PetscInt_FMT " != %" PetscInt_FMT " NumComp2", Nc1, Nc2);
24792df84da0SMatthew G. Knepley 
24802df84da0SMatthew G. Knepley   dim   = dim1 + dim2;
24812df84da0SMatthew G. Knepley   Nc    = Nc1;
24822df84da0SMatthew G. Knepley   Np    = Np1 * Np2;
24832df84da0SMatthew G. Knepley   order = order1;
24849566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q));
24859566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetOrder(*q, order));
24869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Np * dim, &x));
24879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Np, &w));
24882df84da0SMatthew G. Knepley   for (qa = 0, qc = 0; qa < Np1; ++qa) {
24892df84da0SMatthew G. Knepley     for (qb = 0; qb < Np2; ++qb, ++qc) {
2490ad540459SPierre Jolivet       for (d1 = 0, d = 0; d1 < dim1; ++d1, ++d) x[qc * dim + d] = x1[qa * dim1 + d1];
2491ad540459SPierre Jolivet       for (d2 = 0; d2 < dim2; ++d2, ++d) x[qc * dim + d] = x2[qb * dim2 + d2];
24922df84da0SMatthew G. Knepley       w[qc] = w1[qa] * w2[qb];
24932df84da0SMatthew G. Knepley     }
24942df84da0SMatthew G. Knepley   }
24959566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*q, dim, Nc, Np, x, w));
24963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24972df84da0SMatthew G. Knepley }
24982df84da0SMatthew G. Knepley 
2499194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
2500dce8aebaSBarry Smith    A in column-major format
2501dce8aebaSBarry Smith    Ainv in row-major format
2502dce8aebaSBarry Smith    tau has length m
2503dce8aebaSBarry Smith    worksize must be >= max(1,n)
2504194825f6SJed Brown  */
2505d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m, PetscInt mstride, PetscInt n, PetscReal *A_in, PetscReal *Ainv_out, PetscScalar *tau, PetscInt worksize, PetscScalar *work)
2506d71ae5a4SJacob Faibussowitsch {
2507194825f6SJed Brown   PetscBLASInt M, N, K, lda, ldb, ldwork, info;
2508194825f6SJed Brown   PetscScalar *A, *Ainv, *R, *Q, Alpha;
2509194825f6SJed Brown 
2510194825f6SJed Brown   PetscFunctionBegin;
2511194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
2512194825f6SJed Brown   {
2513194825f6SJed Brown     PetscInt i, j;
25149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(m * n, &A, m * n, &Ainv));
2515194825f6SJed Brown     for (j = 0; j < n; j++) {
2516194825f6SJed Brown       for (i = 0; i < m; i++) A[i + m * j] = A_in[i + mstride * j];
2517194825f6SJed Brown     }
2518194825f6SJed Brown     mstride = m;
2519194825f6SJed Brown   }
2520194825f6SJed Brown #else
2521194825f6SJed Brown   A = A_in;
2522194825f6SJed Brown   Ainv = Ainv_out;
2523194825f6SJed Brown #endif
2524194825f6SJed Brown 
25259566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(m, &M));
25269566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &N));
25279566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(mstride, &lda));
25289566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(worksize, &ldwork));
25299566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
2530792fecdfSBarry Smith   PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&M, &N, A, &lda, tau, work, &ldwork, &info));
25319566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPop());
253228b400f6SJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "xGEQRF error");
2533194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
2534194825f6SJed Brown 
2535194825f6SJed Brown   /* Extract an explicit representation of Q */
2536194825f6SJed Brown   Q = Ainv;
25379566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(Q, A, mstride * n));
2538194825f6SJed Brown   K = N; /* full rank */
2539792fecdfSBarry Smith   PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&M, &N, &K, Q, &lda, tau, work, &ldwork, &info));
254028b400f6SJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "xORGQR/xUNGQR error");
2541194825f6SJed Brown 
2542194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
2543194825f6SJed Brown   Alpha = 1.0;
2544194825f6SJed Brown   ldb   = lda;
2545792fecdfSBarry Smith   PetscCallBLAS("BLAStrsm", BLAStrsm_("Right", "Upper", "ConjugateTranspose", "NotUnitTriangular", &M, &N, &Alpha, R, &lda, Q, &ldb));
2546194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
2547194825f6SJed Brown 
2548194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
2549194825f6SJed Brown   {
2550194825f6SJed Brown     PetscInt i;
2551194825f6SJed Brown     for (i = 0; i < m * n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
25529566063dSJacob Faibussowitsch     PetscCall(PetscFree2(A, Ainv));
2553194825f6SJed Brown   }
2554194825f6SJed Brown #endif
25553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2556194825f6SJed Brown }
2557194825f6SJed Brown 
2558194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
2559d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval, const PetscReal *x, PetscInt ndegree, const PetscInt *degrees, PetscBool Transpose, PetscReal *B)
2560d71ae5a4SJacob Faibussowitsch {
2561194825f6SJed Brown   PetscReal *Bv;
2562194825f6SJed Brown   PetscInt   i, j;
2563194825f6SJed Brown 
2564194825f6SJed Brown   PetscFunctionBegin;
25659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((ninterval + 1) * ndegree, &Bv));
2566194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
25679566063dSJacob Faibussowitsch   PetscCall(PetscDTLegendreEval(ninterval + 1, x, ndegree, degrees, Bv, NULL, NULL));
2568194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
2569194825f6SJed Brown   for (i = 0; i < ninterval; i++) {
2570194825f6SJed Brown     for (j = 0; j < ndegree; j++) {
2571194825f6SJed Brown       if (Transpose) B[i + ninterval * j] = Bv[(i + 1) * ndegree + j] - Bv[i * ndegree + j];
2572194825f6SJed Brown       else B[i * ndegree + j] = Bv[(i + 1) * ndegree + j] - Bv[i * ndegree + j];
2573194825f6SJed Brown     }
2574194825f6SJed Brown   }
25759566063dSJacob Faibussowitsch   PetscCall(PetscFree(Bv));
25763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2577194825f6SJed Brown }
2578194825f6SJed Brown 
2579194825f6SJed Brown /*@
2580194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
2581194825f6SJed Brown 
2582194825f6SJed Brown    Not Collective
2583194825f6SJed Brown 
25844165533cSJose E. Roman    Input Parameters:
2585194825f6SJed Brown +  degree - degree of reconstruction polynomial
2586194825f6SJed Brown .  nsource - number of source intervals
2587194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
2588194825f6SJed Brown .  ntarget - number of target intervals
2589194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
2590194825f6SJed Brown 
25914165533cSJose E. Roman    Output Parameter:
2592194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
2593194825f6SJed Brown 
2594194825f6SJed Brown    Level: advanced
2595194825f6SJed Brown 
2596db781477SPatrick Sanan .seealso: `PetscDTLegendreEval()`
2597194825f6SJed Brown @*/
2598d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTReconstructPoly(PetscInt degree, PetscInt nsource, const PetscReal *sourcex, PetscInt ntarget, const PetscReal *targetx, PetscReal *R)
2599d71ae5a4SJacob Faibussowitsch {
2600194825f6SJed Brown   PetscInt     i, j, k, *bdegrees, worksize;
2601194825f6SJed Brown   PetscReal    xmin, xmax, center, hscale, *sourcey, *targety, *Bsource, *Bsinv, *Btarget;
2602194825f6SJed Brown   PetscScalar *tau, *work;
2603194825f6SJed Brown 
2604194825f6SJed Brown   PetscFunctionBegin;
2605194825f6SJed Brown   PetscValidRealPointer(sourcex, 3);
2606194825f6SJed Brown   PetscValidRealPointer(targetx, 5);
2607194825f6SJed Brown   PetscValidRealPointer(R, 6);
260863a3b9bcSJacob Faibussowitsch   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);
260976bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
2610ad540459SPierre Jolivet     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]);
2611ad540459SPierre Jolivet     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]);
261276bd3646SJed Brown   }
2613194825f6SJed Brown   xmin     = PetscMin(sourcex[0], targetx[0]);
2614194825f6SJed Brown   xmax     = PetscMax(sourcex[nsource], targetx[ntarget]);
2615194825f6SJed Brown   center   = (xmin + xmax) / 2;
2616194825f6SJed Brown   hscale   = (xmax - xmin) / 2;
2617194825f6SJed Brown   worksize = nsource;
26189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(degree + 1, &bdegrees, nsource + 1, &sourcey, nsource * (degree + 1), &Bsource, worksize, &work));
26199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(nsource, &tau, nsource * (degree + 1), &Bsinv, ntarget + 1, &targety, ntarget * (degree + 1), &Btarget));
2620194825f6SJed Brown   for (i = 0; i <= nsource; i++) sourcey[i] = (sourcex[i] - center) / hscale;
2621194825f6SJed Brown   for (i = 0; i <= degree; i++) bdegrees[i] = i + 1;
26229566063dSJacob Faibussowitsch   PetscCall(PetscDTLegendreIntegrate(nsource, sourcey, degree + 1, bdegrees, PETSC_TRUE, Bsource));
26239566063dSJacob Faibussowitsch   PetscCall(PetscDTPseudoInverseQR(nsource, nsource, degree + 1, Bsource, Bsinv, tau, nsource, work));
2624194825f6SJed Brown   for (i = 0; i <= ntarget; i++) targety[i] = (targetx[i] - center) / hscale;
26259566063dSJacob Faibussowitsch   PetscCall(PetscDTLegendreIntegrate(ntarget, targety, degree + 1, bdegrees, PETSC_FALSE, Btarget));
2626194825f6SJed Brown   for (i = 0; i < ntarget; i++) {
2627194825f6SJed Brown     PetscReal rowsum = 0;
2628194825f6SJed Brown     for (j = 0; j < nsource; j++) {
2629194825f6SJed Brown       PetscReal sum = 0;
2630ad540459SPierre Jolivet       for (k = 0; k < degree + 1; k++) sum += Btarget[i * (degree + 1) + k] * Bsinv[k * nsource + j];
2631194825f6SJed Brown       R[i * nsource + j] = sum;
2632194825f6SJed Brown       rowsum += sum;
2633194825f6SJed Brown     }
2634194825f6SJed Brown     for (j = 0; j < nsource; j++) R[i * nsource + j] /= rowsum; /* normalize each row */
2635194825f6SJed Brown   }
26369566063dSJacob Faibussowitsch   PetscCall(PetscFree4(bdegrees, sourcey, Bsource, work));
26379566063dSJacob Faibussowitsch   PetscCall(PetscFree4(tau, Bsinv, targety, Btarget));
26383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2639194825f6SJed Brown }
2640916e780bShannah_mairs 
2641916e780bShannah_mairs /*@C
2642916e780bShannah_mairs    PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points
2643916e780bShannah_mairs 
2644916e780bShannah_mairs    Not Collective
2645916e780bShannah_mairs 
2646d8d19677SJose E. Roman    Input Parameters:
2647916e780bShannah_mairs +  n - the number of GLL nodes
2648916e780bShannah_mairs .  nodes - the GLL nodes
2649916e780bShannah_mairs .  weights - the GLL weights
2650f0fc11ceSJed Brown -  f - the function values at the nodes
2651916e780bShannah_mairs 
2652916e780bShannah_mairs    Output Parameter:
2653916e780bShannah_mairs .  in - the value of the integral
2654916e780bShannah_mairs 
2655916e780bShannah_mairs    Level: beginner
2656916e780bShannah_mairs 
2657db781477SPatrick Sanan .seealso: `PetscDTGaussLobattoLegendreQuadrature()`
2658916e780bShannah_mairs @*/
2659d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n, PetscReal *nodes, PetscReal *weights, const PetscReal *f, PetscReal *in)
2660d71ae5a4SJacob Faibussowitsch {
2661916e780bShannah_mairs   PetscInt i;
2662916e780bShannah_mairs 
2663916e780bShannah_mairs   PetscFunctionBegin;
2664916e780bShannah_mairs   *in = 0.;
2665ad540459SPierre Jolivet   for (i = 0; i < n; i++) *in += f[i] * f[i] * weights[i];
26663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2667916e780bShannah_mairs }
2668916e780bShannah_mairs 
2669916e780bShannah_mairs /*@C
2670916e780bShannah_mairs    PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element
2671916e780bShannah_mairs 
2672916e780bShannah_mairs    Not Collective
2673916e780bShannah_mairs 
2674d8d19677SJose E. Roman    Input Parameters:
2675916e780bShannah_mairs +  n - the number of GLL nodes
2676916e780bShannah_mairs .  nodes - the GLL nodes
2677f0fc11ceSJed Brown -  weights - the GLL weights
2678916e780bShannah_mairs 
2679916e780bShannah_mairs    Output Parameter:
2680916e780bShannah_mairs .  A - the stiffness element
2681916e780bShannah_mairs 
2682916e780bShannah_mairs    Level: beginner
2683916e780bShannah_mairs 
2684916e780bShannah_mairs    Notes:
2685dce8aebaSBarry Smith    Destroy this with `PetscGaussLobattoLegendreElementLaplacianDestroy()`
2686916e780bShannah_mairs 
2687916e780bShannah_mairs    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)
2688916e780bShannah_mairs 
2689db781477SPatrick Sanan .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianDestroy()`
2690916e780bShannah_mairs @*/
2691d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA)
2692d71ae5a4SJacob Faibussowitsch {
2693916e780bShannah_mairs   PetscReal      **A;
2694916e780bShannah_mairs   const PetscReal *gllnodes = nodes;
2695916e780bShannah_mairs   const PetscInt   p        = n - 1;
2696916e780bShannah_mairs   PetscReal        z0, z1, z2 = -1, x, Lpj, Lpr;
2697916e780bShannah_mairs   PetscInt         i, j, nn, r;
2698916e780bShannah_mairs 
2699916e780bShannah_mairs   PetscFunctionBegin;
27009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &A));
27019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n * n, &A[0]));
2702916e780bShannah_mairs   for (i = 1; i < n; i++) A[i] = A[i - 1] + n;
2703916e780bShannah_mairs 
2704916e780bShannah_mairs   for (j = 1; j < p; j++) {
2705916e780bShannah_mairs     x  = gllnodes[j];
2706916e780bShannah_mairs     z0 = 1.;
2707916e780bShannah_mairs     z1 = x;
2708916e780bShannah_mairs     for (nn = 1; nn < p; nn++) {
2709916e780bShannah_mairs       z2 = x * z1 * (2. * ((PetscReal)nn) + 1.) / (((PetscReal)nn) + 1.) - z0 * (((PetscReal)nn) / (((PetscReal)nn) + 1.));
2710916e780bShannah_mairs       z0 = z1;
2711916e780bShannah_mairs       z1 = z2;
2712916e780bShannah_mairs     }
2713916e780bShannah_mairs     Lpj = z2;
2714916e780bShannah_mairs     for (r = 1; r < p; r++) {
2715916e780bShannah_mairs       if (r == j) {
2716916e780bShannah_mairs         A[j][j] = 2. / (3. * (1. - gllnodes[j] * gllnodes[j]) * Lpj * Lpj);
2717916e780bShannah_mairs       } else {
2718916e780bShannah_mairs         x  = gllnodes[r];
2719916e780bShannah_mairs         z0 = 1.;
2720916e780bShannah_mairs         z1 = x;
2721916e780bShannah_mairs         for (nn = 1; nn < p; nn++) {
2722916e780bShannah_mairs           z2 = x * z1 * (2. * ((PetscReal)nn) + 1.) / (((PetscReal)nn) + 1.) - z0 * (((PetscReal)nn) / (((PetscReal)nn) + 1.));
2723916e780bShannah_mairs           z0 = z1;
2724916e780bShannah_mairs           z1 = z2;
2725916e780bShannah_mairs         }
2726916e780bShannah_mairs         Lpr     = z2;
2727916e780bShannah_mairs         A[r][j] = 4. / (((PetscReal)p) * (((PetscReal)p) + 1.) * Lpj * Lpr * (gllnodes[j] - gllnodes[r]) * (gllnodes[j] - gllnodes[r]));
2728916e780bShannah_mairs       }
2729916e780bShannah_mairs     }
2730916e780bShannah_mairs   }
2731916e780bShannah_mairs   for (j = 1; j < p + 1; j++) {
2732916e780bShannah_mairs     x  = gllnodes[j];
2733916e780bShannah_mairs     z0 = 1.;
2734916e780bShannah_mairs     z1 = x;
2735916e780bShannah_mairs     for (nn = 1; nn < p; nn++) {
2736916e780bShannah_mairs       z2 = x * z1 * (2. * ((PetscReal)nn) + 1.) / (((PetscReal)nn) + 1.) - z0 * (((PetscReal)nn) / (((PetscReal)nn) + 1.));
2737916e780bShannah_mairs       z0 = z1;
2738916e780bShannah_mairs       z1 = z2;
2739916e780bShannah_mairs     }
2740916e780bShannah_mairs     Lpj     = z2;
2741916e780bShannah_mairs     A[j][0] = 4. * PetscPowRealInt(-1., p) / (((PetscReal)p) * (((PetscReal)p) + 1.) * Lpj * (1. + gllnodes[j]) * (1. + gllnodes[j]));
2742916e780bShannah_mairs     A[0][j] = A[j][0];
2743916e780bShannah_mairs   }
2744916e780bShannah_mairs   for (j = 0; j < p; j++) {
2745916e780bShannah_mairs     x  = gllnodes[j];
2746916e780bShannah_mairs     z0 = 1.;
2747916e780bShannah_mairs     z1 = x;
2748916e780bShannah_mairs     for (nn = 1; nn < p; nn++) {
2749916e780bShannah_mairs       z2 = x * z1 * (2. * ((PetscReal)nn) + 1.) / (((PetscReal)nn) + 1.) - z0 * (((PetscReal)nn) / (((PetscReal)nn) + 1.));
2750916e780bShannah_mairs       z0 = z1;
2751916e780bShannah_mairs       z1 = z2;
2752916e780bShannah_mairs     }
2753916e780bShannah_mairs     Lpj = z2;
2754916e780bShannah_mairs 
2755916e780bShannah_mairs     A[p][j] = 4. / (((PetscReal)p) * (((PetscReal)p) + 1.) * Lpj * (1. - gllnodes[j]) * (1. - gllnodes[j]));
2756916e780bShannah_mairs     A[j][p] = A[p][j];
2757916e780bShannah_mairs   }
2758916e780bShannah_mairs   A[0][0] = 0.5 + (((PetscReal)p) * (((PetscReal)p) + 1.) - 2.) / 6.;
2759916e780bShannah_mairs   A[p][p] = A[0][0];
2760916e780bShannah_mairs   *AA     = A;
27613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2762916e780bShannah_mairs }
2763916e780bShannah_mairs 
2764916e780bShannah_mairs /*@C
2765dce8aebaSBarry Smith    PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element created with `PetscGaussLobattoLegendreElementLaplacianCreate()`
2766916e780bShannah_mairs 
2767916e780bShannah_mairs    Not Collective
2768916e780bShannah_mairs 
2769d8d19677SJose E. Roman    Input Parameters:
2770916e780bShannah_mairs +  n - the number of GLL nodes
2771916e780bShannah_mairs .  nodes - the GLL nodes
2772916e780bShannah_mairs .  weights - the GLL weightss
2773916e780bShannah_mairs -  A - the stiffness element
2774916e780bShannah_mairs 
2775916e780bShannah_mairs    Level: beginner
2776916e780bShannah_mairs 
2777db781477SPatrick Sanan .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianCreate()`
2778916e780bShannah_mairs @*/
2779d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA)
2780d71ae5a4SJacob Faibussowitsch {
2781916e780bShannah_mairs   PetscFunctionBegin;
27829566063dSJacob Faibussowitsch   PetscCall(PetscFree((*AA)[0]));
27839566063dSJacob Faibussowitsch   PetscCall(PetscFree(*AA));
2784916e780bShannah_mairs   *AA = NULL;
27853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2786916e780bShannah_mairs }
2787916e780bShannah_mairs 
2788916e780bShannah_mairs /*@C
2789916e780bShannah_mairs    PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element
2790916e780bShannah_mairs 
2791916e780bShannah_mairs    Not Collective
2792916e780bShannah_mairs 
2793916e780bShannah_mairs    Input Parameter:
2794916e780bShannah_mairs +  n - the number of GLL nodes
2795916e780bShannah_mairs .  nodes - the GLL nodes
2796916e780bShannah_mairs .  weights - the GLL weights
2797916e780bShannah_mairs 
2798d8d19677SJose E. Roman    Output Parameters:
2799916e780bShannah_mairs .  AA - the stiffness element
280020f4b53cSBarry Smith -  AAT - the transpose of AA (pass in `NULL` if you do not need this array)
2801916e780bShannah_mairs 
2802916e780bShannah_mairs    Level: beginner
2803916e780bShannah_mairs 
2804916e780bShannah_mairs    Notes:
2805dce8aebaSBarry Smith    Destroy this with `PetscGaussLobattoLegendreElementGradientDestroy()`
2806916e780bShannah_mairs 
2807916e780bShannah_mairs    You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented
2808916e780bShannah_mairs 
2809dce8aebaSBarry Smith .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianDestroy()`, `PetscGaussLobattoLegendreElementGradientDestroy()`
2810916e780bShannah_mairs @*/
2811d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA, PetscReal ***AAT)
2812d71ae5a4SJacob Faibussowitsch {
2813916e780bShannah_mairs   PetscReal      **A, **AT = NULL;
2814916e780bShannah_mairs   const PetscReal *gllnodes = nodes;
2815916e780bShannah_mairs   const PetscInt   p        = n - 1;
2816e6a796c3SToby Isaac   PetscReal        Li, Lj, d0;
2817916e780bShannah_mairs   PetscInt         i, j;
2818916e780bShannah_mairs 
2819916e780bShannah_mairs   PetscFunctionBegin;
28209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &A));
28219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n * n, &A[0]));
2822916e780bShannah_mairs   for (i = 1; i < n; i++) A[i] = A[i - 1] + n;
2823916e780bShannah_mairs 
2824916e780bShannah_mairs   if (AAT) {
28259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &AT));
28269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n * n, &AT[0]));
2827916e780bShannah_mairs     for (i = 1; i < n; i++) AT[i] = AT[i - 1] + n;
2828916e780bShannah_mairs   }
2829916e780bShannah_mairs 
2830ad540459SPierre Jolivet   if (n == 1) A[0][0] = 0.;
2831916e780bShannah_mairs   d0 = (PetscReal)p * ((PetscReal)p + 1.) / 4.;
2832916e780bShannah_mairs   for (i = 0; i < n; i++) {
2833916e780bShannah_mairs     for (j = 0; j < n; j++) {
2834916e780bShannah_mairs       A[i][j] = 0.;
28359566063dSJacob Faibussowitsch       PetscCall(PetscDTComputeJacobi(0., 0., p, gllnodes[i], &Li));
28369566063dSJacob Faibussowitsch       PetscCall(PetscDTComputeJacobi(0., 0., p, gllnodes[j], &Lj));
2837916e780bShannah_mairs       if (i != j) A[i][j] = Li / (Lj * (gllnodes[i] - gllnodes[j]));
2838916e780bShannah_mairs       if ((j == i) && (i == 0)) A[i][j] = -d0;
2839916e780bShannah_mairs       if (j == i && i == p) A[i][j] = d0;
2840916e780bShannah_mairs       if (AT) AT[j][i] = A[i][j];
2841916e780bShannah_mairs     }
2842916e780bShannah_mairs   }
2843916e780bShannah_mairs   if (AAT) *AAT = AT;
2844916e780bShannah_mairs   *AA = A;
28453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2846916e780bShannah_mairs }
2847916e780bShannah_mairs 
2848916e780bShannah_mairs /*@C
2849dce8aebaSBarry Smith    PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with `PetscGaussLobattoLegendreElementGradientCreate()`
2850916e780bShannah_mairs 
2851916e780bShannah_mairs    Not Collective
2852916e780bShannah_mairs 
2853d8d19677SJose E. Roman    Input Parameters:
2854916e780bShannah_mairs +  n - the number of GLL nodes
2855916e780bShannah_mairs .  nodes - the GLL nodes
2856916e780bShannah_mairs .  weights - the GLL weights
2857916e780bShannah_mairs .  AA - the stiffness element
2858916e780bShannah_mairs -  AAT - the transpose of the element
2859916e780bShannah_mairs 
2860916e780bShannah_mairs    Level: beginner
2861916e780bShannah_mairs 
2862db781477SPatrick Sanan .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianCreate()`, `PetscGaussLobattoLegendreElementAdvectionCreate()`
2863916e780bShannah_mairs @*/
2864d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA, PetscReal ***AAT)
2865d71ae5a4SJacob Faibussowitsch {
2866916e780bShannah_mairs   PetscFunctionBegin;
28679566063dSJacob Faibussowitsch   PetscCall(PetscFree((*AA)[0]));
28689566063dSJacob Faibussowitsch   PetscCall(PetscFree(*AA));
2869916e780bShannah_mairs   *AA = NULL;
2870916e780bShannah_mairs   if (*AAT) {
28719566063dSJacob Faibussowitsch     PetscCall(PetscFree((*AAT)[0]));
28729566063dSJacob Faibussowitsch     PetscCall(PetscFree(*AAT));
2873916e780bShannah_mairs     *AAT = NULL;
2874916e780bShannah_mairs   }
28753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2876916e780bShannah_mairs }
2877916e780bShannah_mairs 
2878916e780bShannah_mairs /*@C
2879916e780bShannah_mairs    PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element
2880916e780bShannah_mairs 
2881916e780bShannah_mairs    Not Collective
2882916e780bShannah_mairs 
2883d8d19677SJose E. Roman    Input Parameters:
2884916e780bShannah_mairs +  n - the number of GLL nodes
2885916e780bShannah_mairs .  nodes - the GLL nodes
2886f0fc11ceSJed Brown -  weights - the GLL weightss
2887916e780bShannah_mairs 
2888916e780bShannah_mairs    Output Parameter:
2889916e780bShannah_mairs .  AA - the stiffness element
2890916e780bShannah_mairs 
2891916e780bShannah_mairs    Level: beginner
2892916e780bShannah_mairs 
2893916e780bShannah_mairs    Notes:
2894dce8aebaSBarry Smith    Destroy this with `PetscGaussLobattoLegendreElementAdvectionDestroy()`
2895916e780bShannah_mairs 
2896916e780bShannah_mairs    This is the same as the Gradient operator multiplied by the diagonal mass matrix
2897916e780bShannah_mairs 
2898916e780bShannah_mairs    You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented
2899916e780bShannah_mairs 
2900db781477SPatrick Sanan .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementLaplacianCreate()`, `PetscGaussLobattoLegendreElementAdvectionDestroy()`
2901916e780bShannah_mairs @*/
2902d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA)
2903d71ae5a4SJacob Faibussowitsch {
2904916e780bShannah_mairs   PetscReal      **D;
2905916e780bShannah_mairs   const PetscReal *gllweights = weights;
2906916e780bShannah_mairs   const PetscInt   glln       = n;
2907916e780bShannah_mairs   PetscInt         i, j;
2908916e780bShannah_mairs 
2909916e780bShannah_mairs   PetscFunctionBegin;
29109566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementGradientCreate(n, nodes, weights, &D, NULL));
2911916e780bShannah_mairs   for (i = 0; i < glln; i++) {
2912ad540459SPierre Jolivet     for (j = 0; j < glln; j++) D[i][j] = gllweights[i] * D[i][j];
2913916e780bShannah_mairs   }
2914916e780bShannah_mairs   *AA = D;
29153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2916916e780bShannah_mairs }
2917916e780bShannah_mairs 
2918916e780bShannah_mairs /*@C
2919dce8aebaSBarry Smith    PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element created with `PetscGaussLobattoLegendreElementAdvectionCreate()`
2920916e780bShannah_mairs 
2921916e780bShannah_mairs    Not Collective
2922916e780bShannah_mairs 
2923d8d19677SJose E. Roman    Input Parameters:
2924916e780bShannah_mairs +  n - the number of GLL nodes
2925916e780bShannah_mairs .  nodes - the GLL nodes
2926916e780bShannah_mairs .  weights - the GLL weights
2927916e780bShannah_mairs -  A - advection
2928916e780bShannah_mairs 
2929916e780bShannah_mairs    Level: beginner
2930916e780bShannah_mairs 
2931db781477SPatrick Sanan .seealso: `PetscDTGaussLobattoLegendreQuadrature()`, `PetscGaussLobattoLegendreElementAdvectionCreate()`
2932916e780bShannah_mairs @*/
2933d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA)
2934d71ae5a4SJacob Faibussowitsch {
2935916e780bShannah_mairs   PetscFunctionBegin;
29369566063dSJacob Faibussowitsch   PetscCall(PetscFree((*AA)[0]));
29379566063dSJacob Faibussowitsch   PetscCall(PetscFree(*AA));
2938916e780bShannah_mairs   *AA = NULL;
29393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2940916e780bShannah_mairs }
2941916e780bShannah_mairs 
2942d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA)
2943d71ae5a4SJacob Faibussowitsch {
2944916e780bShannah_mairs   PetscReal      **A;
2945916e780bShannah_mairs   const PetscReal *gllweights = weights;
2946916e780bShannah_mairs   const PetscInt   glln       = n;
2947916e780bShannah_mairs   PetscInt         i, j;
2948916e780bShannah_mairs 
2949916e780bShannah_mairs   PetscFunctionBegin;
29509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(glln, &A));
29519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(glln * glln, &A[0]));
2952916e780bShannah_mairs   for (i = 1; i < glln; i++) A[i] = A[i - 1] + glln;
2953ad540459SPierre Jolivet   if (glln == 1) A[0][0] = 0.;
2954916e780bShannah_mairs   for (i = 0; i < glln; i++) {
2955916e780bShannah_mairs     for (j = 0; j < glln; j++) {
2956916e780bShannah_mairs       A[i][j] = 0.;
2957916e780bShannah_mairs       if (j == i) A[i][j] = gllweights[i];
2958916e780bShannah_mairs     }
2959916e780bShannah_mairs   }
2960916e780bShannah_mairs   *AA = A;
29613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2962916e780bShannah_mairs }
2963916e780bShannah_mairs 
2964d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n, PetscReal *nodes, PetscReal *weights, PetscReal ***AA)
2965d71ae5a4SJacob Faibussowitsch {
2966916e780bShannah_mairs   PetscFunctionBegin;
29679566063dSJacob Faibussowitsch   PetscCall(PetscFree((*AA)[0]));
29689566063dSJacob Faibussowitsch   PetscCall(PetscFree(*AA));
2969916e780bShannah_mairs   *AA = NULL;
29703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2971916e780bShannah_mairs }
2972d4afb720SToby Isaac 
2973d4afb720SToby Isaac /*@
2974d4afb720SToby Isaac   PetscDTIndexToBary - convert an index into a barycentric coordinate.
2975d4afb720SToby Isaac 
2976d4afb720SToby Isaac   Input Parameters:
2977d4afb720SToby Isaac + 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)
2978d4afb720SToby Isaac . sum - the value that the sum of the barycentric coordinates (which will be non-negative integers) should sum to
2979d4afb720SToby Isaac - index - the index to convert: should be >= 0 and < Binomial(len - 1 + sum, sum)
2980d4afb720SToby Isaac 
2981d4afb720SToby Isaac   Output Parameter:
2982d4afb720SToby Isaac . coord - will be filled with the barycentric coordinate
2983d4afb720SToby Isaac 
2984d4afb720SToby Isaac   Level: beginner
2985d4afb720SToby Isaac 
2986dce8aebaSBarry Smith   Note:
2987dce8aebaSBarry Smith   The indices map to barycentric coordinates in lexicographic order, where the first index is the
2988d4afb720SToby Isaac   least significant and the last index is the most significant.
2989d4afb720SToby Isaac 
2990db781477SPatrick Sanan .seealso: `PetscDTBaryToIndex()`
2991d4afb720SToby Isaac @*/
2992d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTIndexToBary(PetscInt len, PetscInt sum, PetscInt index, PetscInt coord[])
2993d71ae5a4SJacob Faibussowitsch {
2994d4afb720SToby Isaac   PetscInt c, d, s, total, subtotal, nexttotal;
2995d4afb720SToby Isaac 
2996d4afb720SToby Isaac   PetscFunctionBeginHot;
299708401ef6SPierre Jolivet   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
299808401ef6SPierre Jolivet   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index must be non-negative");
2999d4afb720SToby Isaac   if (!len) {
30003ba16761SJacob Faibussowitsch     if (!sum && !index) PetscFunctionReturn(PETSC_SUCCESS);
3001d4afb720SToby Isaac     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid index or sum for length 0 barycentric coordinate");
3002d4afb720SToby Isaac   }
3003d4afb720SToby Isaac   for (c = 1, total = 1; c <= len; c++) {
3004d4afb720SToby Isaac     /* total is the number of ways to have a tuple of length c with sum */
3005d4afb720SToby Isaac     if (index < total) break;
3006d4afb720SToby Isaac     total = (total * (sum + c)) / c;
3007d4afb720SToby Isaac   }
300808401ef6SPierre Jolivet   PetscCheck(c <= len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index out of range");
3009d4afb720SToby Isaac   for (d = c; d < len; d++) coord[d] = 0;
3010d4afb720SToby Isaac   for (s = 0, subtotal = 1, nexttotal = 1; c > 0;) {
3011d4afb720SToby Isaac     /* subtotal is the number of ways to have a tuple of length c with sum s */
3012d4afb720SToby Isaac     /* nexttotal is the number of ways to have a tuple of length c-1 with sum s */
3013d4afb720SToby Isaac     if ((index + subtotal) >= total) {
3014d4afb720SToby Isaac       coord[--c] = sum - s;
3015d4afb720SToby Isaac       index -= (total - subtotal);
3016d4afb720SToby Isaac       sum       = s;
3017d4afb720SToby Isaac       total     = nexttotal;
3018d4afb720SToby Isaac       subtotal  = 1;
3019d4afb720SToby Isaac       nexttotal = 1;
3020d4afb720SToby Isaac       s         = 0;
3021d4afb720SToby Isaac     } else {
3022d4afb720SToby Isaac       subtotal  = (subtotal * (c + s)) / (s + 1);
3023d4afb720SToby Isaac       nexttotal = (nexttotal * (c - 1 + s)) / (s + 1);
3024d4afb720SToby Isaac       s++;
3025d4afb720SToby Isaac     }
3026d4afb720SToby Isaac   }
30273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3028d4afb720SToby Isaac }
3029d4afb720SToby Isaac 
3030d4afb720SToby Isaac /*@
3031d4afb720SToby Isaac   PetscDTBaryToIndex - convert a barycentric coordinate to an index
3032d4afb720SToby Isaac 
3033d4afb720SToby Isaac   Input Parameters:
3034d4afb720SToby Isaac + 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)
3035d4afb720SToby Isaac . sum - the value that the sum of the barycentric coordinates (which will be non-negative integers) should sum to
3036d4afb720SToby Isaac - coord - a barycentric coordinate with the given length and sum
3037d4afb720SToby Isaac 
3038d4afb720SToby Isaac   Output Parameter:
3039d4afb720SToby Isaac . index - the unique index for the coordinate, >= 0 and < Binomial(len - 1 + sum, sum)
3040d4afb720SToby Isaac 
3041d4afb720SToby Isaac   Level: beginner
3042d4afb720SToby Isaac 
3043dce8aebaSBarry Smith   Note:
3044dce8aebaSBarry Smith   The indices map to barycentric coordinates in lexicographic order, where the first index is the
3045d4afb720SToby Isaac   least significant and the last index is the most significant.
3046d4afb720SToby Isaac 
3047db781477SPatrick Sanan .seealso: `PetscDTIndexToBary`
3048d4afb720SToby Isaac @*/
3049d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTBaryToIndex(PetscInt len, PetscInt sum, const PetscInt coord[], PetscInt *index)
3050d71ae5a4SJacob Faibussowitsch {
3051d4afb720SToby Isaac   PetscInt c;
3052d4afb720SToby Isaac   PetscInt i;
3053d4afb720SToby Isaac   PetscInt total;
3054d4afb720SToby Isaac 
3055d4afb720SToby Isaac   PetscFunctionBeginHot;
305608401ef6SPierre Jolivet   PetscCheck(len >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative");
3057d4afb720SToby Isaac   if (!len) {
3058d4afb720SToby Isaac     if (!sum) {
3059d4afb720SToby Isaac       *index = 0;
30603ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
3061d4afb720SToby Isaac     }
3062d4afb720SToby Isaac     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid index or sum for length 0 barycentric coordinate");
3063d4afb720SToby Isaac   }
3064d4afb720SToby Isaac   for (c = 1, total = 1; c < len; c++) total = (total * (sum + c)) / c;
3065d4afb720SToby Isaac   i = total - 1;
3066d4afb720SToby Isaac   c = len - 1;
3067d4afb720SToby Isaac   sum -= coord[c];
3068d4afb720SToby Isaac   while (sum > 0) {
3069d4afb720SToby Isaac     PetscInt subtotal;
3070d4afb720SToby Isaac     PetscInt s;
3071d4afb720SToby Isaac 
3072d4afb720SToby Isaac     for (s = 1, subtotal = 1; s < sum; s++) subtotal = (subtotal * (c + s)) / s;
3073d4afb720SToby Isaac     i -= subtotal;
3074d4afb720SToby Isaac     sum -= coord[--c];
3075d4afb720SToby Isaac   }
3076d4afb720SToby Isaac   *index = i;
30773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3078d4afb720SToby Isaac }
3079*07218a29SMatthew G. Knepley 
3080*07218a29SMatthew G. Knepley /*
3081*07218a29SMatthew G. Knepley 1) For each face type:
3082*07218a29SMatthew G. Knepley      For each transformation from outward to inward normal:
3083*07218a29SMatthew G. Knepley          Compute the permutation of quadrature points:
3084*07218a29SMatthew G. Knepley             Compute the quad point coordinates from each side and compare
3085*07218a29SMatthew G. Knepley */
3086*07218a29SMatthew G. Knepley PetscErrorCode PetscDTComputeFaceQuadPermutation(DMPolytopeType ct, PetscQuadrature quad, PetscInt *Np, IS *perm[])
3087*07218a29SMatthew G. Knepley {
3088*07218a29SMatthew G. Knepley   const PetscReal *xq, *wq;
3089*07218a29SMatthew G. Knepley   PetscInt         dim, qdim, d, Na, o, Nq, q, qp;
3090*07218a29SMatthew G. Knepley 
3091*07218a29SMatthew G. Knepley   PetscFunctionBegin;
3092*07218a29SMatthew G. Knepley   dim = DMPolytopeTypeGetDim(ct);
3093*07218a29SMatthew G. Knepley   Na  = DMPolytopeTypeGetNumArrangments(ct);
3094*07218a29SMatthew G. Knepley   Na /= 2;
3095*07218a29SMatthew G. Knepley   PetscCall(PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &xq, &wq));
3096*07218a29SMatthew G. Knepley   PetscCheck(dim >= 0 && dim < 3, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot compute quadrature permutation for cell type %s", DMPolytopeTypes[ct]);
3097*07218a29SMatthew G. Knepley   PetscCheck(dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Celltype %s dimension %" PetscInt_FMT " != %" PetscInt_FMT " quad dimension", DMPolytopeTypes[ct], dim, qdim);
3098*07218a29SMatthew G. Knepley   *Np = Na;
3099*07218a29SMatthew G. Knepley   PetscCall(PetscMalloc1(Na, perm));
3100*07218a29SMatthew G. Knepley   for (o = -1; o > -(Na + 1); --o) {
3101*07218a29SMatthew G. Knepley     DM        refdm;
3102*07218a29SMatthew G. Knepley     PetscInt *idx;
3103*07218a29SMatthew G. Knepley     PetscReal xi0[3] = {-1., -1., -1.}, v0[3], J[9], detJ, txq[3];
3104*07218a29SMatthew G. Knepley     PetscBool flg;
3105*07218a29SMatthew G. Knepley 
3106*07218a29SMatthew G. Knepley     PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &refdm));
3107*07218a29SMatthew G. Knepley     PetscCall(DMPlexOrientPoint(refdm, 0, o));
3108*07218a29SMatthew G. Knepley     PetscCall(DMPlexComputeCellGeometryFEM(refdm, 0, NULL, v0, J, NULL, &detJ));
3109*07218a29SMatthew G. Knepley     PetscCall(PetscMalloc1(Nq, &idx));
3110*07218a29SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
3111*07218a29SMatthew G. Knepley       CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], txq);
3112*07218a29SMatthew G. Knepley       for (qp = 0; qp < Nq; ++qp) {
3113*07218a29SMatthew G. Knepley         PetscReal diff = 0.;
3114*07218a29SMatthew G. Knepley 
3115*07218a29SMatthew G. Knepley         for (d = 0; d < dim; ++d) diff += PetscAbsReal(txq[d] - xq[qp * dim + d]);
3116*07218a29SMatthew G. Knepley         if (diff < PETSC_SMALL) break;
3117*07218a29SMatthew G. Knepley       }
3118*07218a29SMatthew G. Knepley       PetscCheck(qp < Nq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Transformed quad point %" PetscInt_FMT " does not match another quad point", q);
3119*07218a29SMatthew G. Knepley       idx[q] = qp;
3120*07218a29SMatthew G. Knepley     }
3121*07218a29SMatthew G. Knepley     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nq, idx, PETSC_OWN_POINTER, &(*perm)[-(o + 1)]));
3122*07218a29SMatthew G. Knepley     PetscCall(ISGetInfo((*perm)[-(o + 1)], IS_PERMUTATION, IS_LOCAL, PETSC_TRUE, &flg));
3123*07218a29SMatthew G. Knepley     PetscCall(DMDestroy(&refdm));
3124*07218a29SMatthew G. Knepley     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ordering for orientation %" PetscInt_FMT "was not a permutation", o);
3125*07218a29SMatthew G. Knepley   }
3126*07218a29SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3127*07218a29SMatthew G. Knepley }
3128