xref: /petsc/src/dm/dt/interface/dt.c (revision fbfcfee5110779e3d6a9465ca0a2e0f9a1a6e5e3)
137045ce4SJed Brown /* Discretization tools */
237045ce4SJed Brown 
3a6fc04d9SSatish Balay #include <petscconf.h>
4a6fc04d9SSatish Balay #if defined(PETSC_HAVE_MATHIMF_H)
5a6fc04d9SSatish Balay #include <mathimf.h>           /* this needs to be included before math.h */
6a6fc04d9SSatish Balay #endif
729f144ccSMatthew G. Knepley #ifdef PETSC_HAVE_MPFR
829f144ccSMatthew G. Knepley #include <mpfr.h>
929f144ccSMatthew G. Knepley #endif
10a6fc04d9SSatish Balay 
110c35b76eSJed Brown #include <petscdt.h>            /*I "petscdt.h" I*/
1237045ce4SJed Brown #include <petscblaslapack.h>
13af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
14af0996ceSBarry Smith #include <petsc/private/dtimpl.h>
15665c2dedSJed Brown #include <petscviewer.h>
1659804f93SMatthew G. Knepley #include <petscdmplex.h>
1759804f93SMatthew G. Knepley #include <petscdmshell.h>
1837045ce4SJed Brown 
190bfcf5a5SMatthew G. Knepley static PetscBool GaussCite       = PETSC_FALSE;
200bfcf5a5SMatthew G. Knepley const char       GaussCitation[] = "@article{GolubWelsch1969,\n"
210bfcf5a5SMatthew G. Knepley                                    "  author  = {Golub and Welsch},\n"
220bfcf5a5SMatthew G. Knepley                                    "  title   = {Calculation of Quadrature Rules},\n"
230bfcf5a5SMatthew G. Knepley                                    "  journal = {Math. Comp.},\n"
240bfcf5a5SMatthew G. Knepley                                    "  volume  = {23},\n"
250bfcf5a5SMatthew G. Knepley                                    "  number  = {106},\n"
260bfcf5a5SMatthew G. Knepley                                    "  pages   = {221--230},\n"
270bfcf5a5SMatthew G. Knepley                                    "  year    = {1969}\n}\n";
280bfcf5a5SMatthew G. Knepley 
2940d8ff71SMatthew G. Knepley /*@
3040d8ff71SMatthew G. Knepley   PetscQuadratureCreate - Create a PetscQuadrature object
3140d8ff71SMatthew G. Knepley 
3240d8ff71SMatthew G. Knepley   Collective on MPI_Comm
3340d8ff71SMatthew G. Knepley 
3440d8ff71SMatthew G. Knepley   Input Parameter:
3540d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object
3640d8ff71SMatthew G. Knepley 
3740d8ff71SMatthew G. Knepley   Output Parameter:
3840d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
3940d8ff71SMatthew G. Knepley 
4040d8ff71SMatthew G. Knepley   Level: beginner
4140d8ff71SMatthew G. Knepley 
4240d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, create
4340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
4440d8ff71SMatthew G. Knepley @*/
4521454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
4621454ff5SMatthew G. Knepley {
4721454ff5SMatthew G. Knepley   PetscErrorCode ierr;
4821454ff5SMatthew G. Knepley 
4921454ff5SMatthew G. Knepley   PetscFunctionBegin;
5021454ff5SMatthew G. Knepley   PetscValidPointer(q, 2);
51623436dcSMatthew G. Knepley   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
5273107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*q,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
5321454ff5SMatthew G. Knepley   (*q)->dim       = -1;
54bcede257SMatthew G. Knepley   (*q)->order     = -1;
5521454ff5SMatthew G. Knepley   (*q)->numPoints = 0;
5621454ff5SMatthew G. Knepley   (*q)->points    = NULL;
5721454ff5SMatthew G. Knepley   (*q)->weights   = NULL;
5821454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
5921454ff5SMatthew G. Knepley }
6021454ff5SMatthew G. Knepley 
61c9638911SMatthew G. Knepley /*@
62c9638911SMatthew G. Knepley   PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object
63c9638911SMatthew G. Knepley 
64c9638911SMatthew G. Knepley   Collective on PetscQuadrature
65c9638911SMatthew G. Knepley 
66c9638911SMatthew G. Knepley   Input Parameter:
67c9638911SMatthew G. Knepley . q  - The PetscQuadrature object
68c9638911SMatthew G. Knepley 
69c9638911SMatthew G. Knepley   Output Parameter:
70c9638911SMatthew G. Knepley . r  - The new PetscQuadrature object
71c9638911SMatthew G. Knepley 
72c9638911SMatthew G. Knepley   Level: beginner
73c9638911SMatthew G. Knepley 
74c9638911SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, clone
75c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData()
76c9638911SMatthew G. Knepley @*/
77c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
78c9638911SMatthew G. Knepley {
79c9638911SMatthew G. Knepley   PetscInt         order, dim, Nq;
80c9638911SMatthew G. Knepley   const PetscReal *points, *weights;
81c9638911SMatthew G. Knepley   PetscReal       *p, *w;
82c9638911SMatthew G. Knepley   PetscErrorCode   ierr;
83c9638911SMatthew G. Knepley 
84c9638911SMatthew G. Knepley   PetscFunctionBegin;
85c9638911SMatthew G. Knepley   PetscValidPointer(q, 2);
86c9638911SMatthew G. Knepley   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr);
87c9638911SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
88c9638911SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr);
89c9638911SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &Nq, &points, &weights);CHKERRQ(ierr);
90c9638911SMatthew G. Knepley   ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr);
91c9638911SMatthew G. Knepley   ierr = PetscMalloc1(Nq, &w);CHKERRQ(ierr);
92c9638911SMatthew G. Knepley   ierr = PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));CHKERRQ(ierr);
93c9638911SMatthew G. Knepley   ierr = PetscMemcpy(w, weights, Nq * sizeof(PetscReal));CHKERRQ(ierr);
94c9638911SMatthew G. Knepley   ierr = PetscQuadratureSetData(*r, dim, Nq, p, w);CHKERRQ(ierr);
95c9638911SMatthew G. Knepley   PetscFunctionReturn(0);
96c9638911SMatthew G. Knepley }
97c9638911SMatthew G. Knepley 
9840d8ff71SMatthew G. Knepley /*@
9940d8ff71SMatthew G. Knepley   PetscQuadratureDestroy - Destroys a PetscQuadrature object
10040d8ff71SMatthew G. Knepley 
10140d8ff71SMatthew G. Knepley   Collective on PetscQuadrature
10240d8ff71SMatthew G. Knepley 
10340d8ff71SMatthew G. Knepley   Input Parameter:
10440d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
10540d8ff71SMatthew G. Knepley 
10640d8ff71SMatthew G. Knepley   Level: beginner
10740d8ff71SMatthew G. Knepley 
10840d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, destroy
10940d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
11040d8ff71SMatthew G. Knepley @*/
111bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
112bfa639d9SMatthew G. Knepley {
113bfa639d9SMatthew G. Knepley   PetscErrorCode ierr;
114bfa639d9SMatthew G. Knepley 
115bfa639d9SMatthew G. Knepley   PetscFunctionBegin;
11621454ff5SMatthew G. Knepley   if (!*q) PetscFunctionReturn(0);
11721454ff5SMatthew G. Knepley   PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1);
11821454ff5SMatthew G. Knepley   if (--((PetscObject)(*q))->refct > 0) {
11921454ff5SMatthew G. Knepley     *q = NULL;
12021454ff5SMatthew G. Knepley     PetscFunctionReturn(0);
12121454ff5SMatthew G. Knepley   }
12221454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
12321454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
12421454ff5SMatthew G. Knepley   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
12521454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
12621454ff5SMatthew G. Knepley }
12721454ff5SMatthew G. Knepley 
128bcede257SMatthew G. Knepley /*@
129bcede257SMatthew G. Knepley   PetscQuadratureGetOrder - Return the quadrature information
130bcede257SMatthew G. Knepley 
131bcede257SMatthew G. Knepley   Not collective
132bcede257SMatthew G. Knepley 
133bcede257SMatthew G. Knepley   Input Parameter:
134bcede257SMatthew G. Knepley . q - The PetscQuadrature object
135bcede257SMatthew G. Knepley 
136bcede257SMatthew G. Knepley   Output Parameter:
137bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
138bcede257SMatthew G. Knepley 
139bcede257SMatthew G. Knepley   Output Parameter:
140bcede257SMatthew G. Knepley 
141bcede257SMatthew G. Knepley   Level: intermediate
142bcede257SMatthew G. Knepley 
143bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
144bcede257SMatthew G. Knepley @*/
145bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
146bcede257SMatthew G. Knepley {
147bcede257SMatthew G. Knepley   PetscFunctionBegin;
148bcede257SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
149bcede257SMatthew G. Knepley   PetscValidPointer(order, 2);
150bcede257SMatthew G. Knepley   *order = q->order;
151bcede257SMatthew G. Knepley   PetscFunctionReturn(0);
152bcede257SMatthew G. Knepley }
153bcede257SMatthew G. Knepley 
154bcede257SMatthew G. Knepley /*@
155bcede257SMatthew G. Knepley   PetscQuadratureSetOrder - Return the quadrature information
156bcede257SMatthew G. Knepley 
157bcede257SMatthew G. Knepley   Not collective
158bcede257SMatthew G. Knepley 
159bcede257SMatthew G. Knepley   Input Parameters:
160bcede257SMatthew G. Knepley + q - The PetscQuadrature object
161bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
162bcede257SMatthew G. Knepley 
163bcede257SMatthew G. Knepley   Level: intermediate
164bcede257SMatthew G. Knepley 
165bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
166bcede257SMatthew G. Knepley @*/
167bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
168bcede257SMatthew G. Knepley {
169bcede257SMatthew G. Knepley   PetscFunctionBegin;
170bcede257SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
171bcede257SMatthew G. Knepley   q->order = order;
172bcede257SMatthew G. Knepley   PetscFunctionReturn(0);
173bcede257SMatthew G. Knepley }
174bcede257SMatthew G. Knepley 
17540d8ff71SMatthew G. Knepley /*@C
17640d8ff71SMatthew G. Knepley   PetscQuadratureGetData - Returns the data defining the quadrature
17740d8ff71SMatthew G. Knepley 
17840d8ff71SMatthew G. Knepley   Not collective
17940d8ff71SMatthew G. Knepley 
18040d8ff71SMatthew G. Knepley   Input Parameter:
18140d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
18240d8ff71SMatthew G. Knepley 
18340d8ff71SMatthew G. Knepley   Output Parameters:
18440d8ff71SMatthew G. Knepley + dim - The spatial dimension
18540d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
18640d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
18740d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
18840d8ff71SMatthew G. Knepley 
18940d8ff71SMatthew G. Knepley   Level: intermediate
19040d8ff71SMatthew G. Knepley 
19140d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature
19240d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
19340d8ff71SMatthew G. Knepley @*/
19421454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
19521454ff5SMatthew G. Knepley {
19621454ff5SMatthew G. Knepley   PetscFunctionBegin;
19721454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
19821454ff5SMatthew G. Knepley   if (dim) {
19921454ff5SMatthew G. Knepley     PetscValidPointer(dim, 2);
20021454ff5SMatthew G. Knepley     *dim = q->dim;
20121454ff5SMatthew G. Knepley   }
20221454ff5SMatthew G. Knepley   if (npoints) {
20321454ff5SMatthew G. Knepley     PetscValidPointer(npoints, 3);
20421454ff5SMatthew G. Knepley     *npoints = q->numPoints;
20521454ff5SMatthew G. Knepley   }
20621454ff5SMatthew G. Knepley   if (points) {
20721454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
20821454ff5SMatthew G. Knepley     *points = q->points;
20921454ff5SMatthew G. Knepley   }
21021454ff5SMatthew G. Knepley   if (weights) {
21121454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
21221454ff5SMatthew G. Knepley     *weights = q->weights;
21321454ff5SMatthew G. Knepley   }
21421454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
21521454ff5SMatthew G. Knepley }
21621454ff5SMatthew G. Knepley 
21740d8ff71SMatthew G. Knepley /*@C
21840d8ff71SMatthew G. Knepley   PetscQuadratureSetData - Sets the data defining the quadrature
21940d8ff71SMatthew G. Knepley 
22040d8ff71SMatthew G. Knepley   Not collective
22140d8ff71SMatthew G. Knepley 
22240d8ff71SMatthew G. Knepley   Input Parameters:
22340d8ff71SMatthew G. Knepley + q  - The PetscQuadrature object
22440d8ff71SMatthew G. Knepley . dim - The spatial dimension
22540d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
22640d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
22740d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
22840d8ff71SMatthew G. Knepley 
22940d8ff71SMatthew G. Knepley   Level: intermediate
23040d8ff71SMatthew G. Knepley 
23140d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature
23240d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
23340d8ff71SMatthew G. Knepley @*/
23421454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
23521454ff5SMatthew G. Knepley {
23621454ff5SMatthew G. Knepley   PetscFunctionBegin;
23721454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
23821454ff5SMatthew G. Knepley   if (dim >= 0)     q->dim       = dim;
23921454ff5SMatthew G. Knepley   if (npoints >= 0) q->numPoints = npoints;
24021454ff5SMatthew G. Knepley   if (points) {
24121454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
24221454ff5SMatthew G. Knepley     q->points = points;
24321454ff5SMatthew G. Knepley   }
24421454ff5SMatthew G. Knepley   if (weights) {
24521454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
24621454ff5SMatthew G. Knepley     q->weights = weights;
24721454ff5SMatthew G. Knepley   }
248f9fd7fdbSMatthew G. Knepley   PetscFunctionReturn(0);
249f9fd7fdbSMatthew G. Knepley }
250f9fd7fdbSMatthew G. Knepley 
25140d8ff71SMatthew G. Knepley /*@C
25240d8ff71SMatthew G. Knepley   PetscQuadratureView - Views a PetscQuadrature object
25340d8ff71SMatthew G. Knepley 
25440d8ff71SMatthew G. Knepley   Collective on PetscQuadrature
25540d8ff71SMatthew G. Knepley 
25640d8ff71SMatthew G. Knepley   Input Parameters:
25740d8ff71SMatthew G. Knepley + q  - The PetscQuadrature object
25840d8ff71SMatthew G. Knepley - viewer - The PetscViewer object
25940d8ff71SMatthew G. Knepley 
26040d8ff71SMatthew G. Knepley   Level: beginner
26140d8ff71SMatthew G. Knepley 
26240d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, view
26340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
26440d8ff71SMatthew G. Knepley @*/
265f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
266f9fd7fdbSMatthew G. Knepley {
267f9fd7fdbSMatthew G. Knepley   PetscInt       q, d;
268f9fd7fdbSMatthew G. Knepley   PetscErrorCode ierr;
269f9fd7fdbSMatthew G. Knepley 
270f9fd7fdbSMatthew G. Knepley   PetscFunctionBegin;
27198c3331eSBarry Smith   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr);
27221454ff5SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n  (", quad->numPoints);CHKERRQ(ierr);
27321454ff5SMatthew G. Knepley   for (q = 0; q < quad->numPoints; ++q) {
27421454ff5SMatthew G. Knepley     for (d = 0; d < quad->dim; ++d) {
275f9fd7fdbSMatthew G. Knepley       if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);
276ab15ae43SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
277f9fd7fdbSMatthew G. Knepley     }
278ab15ae43SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr);
279f9fd7fdbSMatthew G. Knepley   }
280bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
281bfa639d9SMatthew G. Knepley }
282bfa639d9SMatthew G. Knepley 
28389710940SMatthew G. Knepley /*@C
28489710940SMatthew G. Knepley   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
28589710940SMatthew G. Knepley 
28689710940SMatthew G. Knepley   Not collective
28789710940SMatthew G. Knepley 
28889710940SMatthew G. Knepley   Input Parameter:
28989710940SMatthew G. Knepley + q - The original PetscQuadrature
29089710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into
29189710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement
29289710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement
29389710940SMatthew G. Knepley 
29489710940SMatthew G. Knepley   Output Parameters:
29589710940SMatthew G. Knepley . dim - The dimension
29689710940SMatthew G. Knepley 
29789710940SMatthew G. Knepley   Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
29889710940SMatthew G. Knepley 
29989710940SMatthew G. Knepley   Level: intermediate
30089710940SMatthew G. Knepley 
30189710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
30289710940SMatthew G. Knepley @*/
30389710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
30489710940SMatthew G. Knepley {
30589710940SMatthew G. Knepley   const PetscReal *points,    *weights;
30689710940SMatthew G. Knepley   PetscReal       *pointsRef, *weightsRef;
30789710940SMatthew G. Knepley   PetscInt         dim, order, npoints, npointsRef, c, p, d, e;
30889710940SMatthew G. Knepley   PetscErrorCode   ierr;
30989710940SMatthew G. Knepley 
31089710940SMatthew G. Knepley   PetscFunctionBegin;
31189710940SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
31289710940SMatthew G. Knepley   PetscValidPointer(v0, 3);
31389710940SMatthew G. Knepley   PetscValidPointer(jac, 4);
31489710940SMatthew G. Knepley   PetscValidPointer(qref, 5);
31589710940SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr);
31689710940SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
31789710940SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &npoints, &points, &weights);CHKERRQ(ierr);
31889710940SMatthew G. Knepley   npointsRef = npoints*numSubelements;
31989710940SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr);
32089710940SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef,&weightsRef);CHKERRQ(ierr);
32189710940SMatthew G. Knepley   for (c = 0; c < numSubelements; ++c) {
32289710940SMatthew G. Knepley     for (p = 0; p < npoints; ++p) {
32389710940SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
32489710940SMatthew G. Knepley         pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
32589710940SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
32689710940SMatthew G. Knepley           pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
32789710940SMatthew G. Knepley         }
32889710940SMatthew G. Knepley       }
32989710940SMatthew G. Knepley       /* Could also use detJ here */
33089710940SMatthew G. Knepley       weightsRef[c*npoints+p] = weights[p]/numSubelements;
33189710940SMatthew G. Knepley     }
33289710940SMatthew G. Knepley   }
33389710940SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr);
33489710940SMatthew G. Knepley   ierr = PetscQuadratureSetData(*qref, dim, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr);
33589710940SMatthew G. Knepley   PetscFunctionReturn(0);
33689710940SMatthew G. Knepley }
33789710940SMatthew G. Knepley 
33837045ce4SJed Brown /*@
33937045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
34037045ce4SJed Brown 
34137045ce4SJed Brown    Not Collective
34237045ce4SJed Brown 
34337045ce4SJed Brown    Input Arguments:
34437045ce4SJed Brown +  npoints - number of spatial points to evaluate at
34537045ce4SJed Brown .  points - array of locations to evaluate at
34637045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
34737045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
34837045ce4SJed Brown 
34937045ce4SJed Brown    Output Arguments:
3500298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
3510298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
3520298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
35337045ce4SJed Brown 
35437045ce4SJed Brown    Level: intermediate
35537045ce4SJed Brown 
35637045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
35737045ce4SJed Brown @*/
35837045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
35937045ce4SJed Brown {
36037045ce4SJed Brown   PetscInt i,maxdegree;
36137045ce4SJed Brown 
36237045ce4SJed Brown   PetscFunctionBegin;
36337045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
36437045ce4SJed Brown   maxdegree = degrees[ndegree-1];
36537045ce4SJed Brown   for (i=0; i<npoints; i++) {
36637045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
36737045ce4SJed Brown     PetscInt  j,k;
36837045ce4SJed Brown     x    = points[i];
36937045ce4SJed Brown     pm2  = 0;
37037045ce4SJed Brown     pm1  = 1;
37137045ce4SJed Brown     pd2  = 0;
37237045ce4SJed Brown     pd1  = 0;
37337045ce4SJed Brown     pdd2 = 0;
37437045ce4SJed Brown     pdd1 = 0;
37537045ce4SJed Brown     k    = 0;
37637045ce4SJed Brown     if (degrees[k] == 0) {
37737045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
37837045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
37937045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
38037045ce4SJed Brown       k++;
38137045ce4SJed Brown     }
38237045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
38337045ce4SJed Brown       PetscReal p,d,dd;
38437045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
38537045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
38637045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
38737045ce4SJed Brown       pm2  = pm1;
38837045ce4SJed Brown       pm1  = p;
38937045ce4SJed Brown       pd2  = pd1;
39037045ce4SJed Brown       pd1  = d;
39137045ce4SJed Brown       pdd2 = pdd1;
39237045ce4SJed Brown       pdd1 = dd;
39337045ce4SJed Brown       if (degrees[k] == j) {
39437045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
39537045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
39637045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
39737045ce4SJed Brown       }
39837045ce4SJed Brown     }
39937045ce4SJed Brown   }
40037045ce4SJed Brown   PetscFunctionReturn(0);
40137045ce4SJed Brown }
40237045ce4SJed Brown 
40337045ce4SJed Brown /*@
40437045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
40537045ce4SJed Brown 
40637045ce4SJed Brown    Not Collective
40737045ce4SJed Brown 
40837045ce4SJed Brown    Input Arguments:
40937045ce4SJed Brown +  npoints - number of points
41037045ce4SJed Brown .  a - left end of interval (often-1)
41137045ce4SJed Brown -  b - right end of interval (often +1)
41237045ce4SJed Brown 
41337045ce4SJed Brown    Output Arguments:
41437045ce4SJed Brown +  x - quadrature points
41537045ce4SJed Brown -  w - quadrature weights
41637045ce4SJed Brown 
41737045ce4SJed Brown    Level: intermediate
41837045ce4SJed Brown 
41937045ce4SJed Brown    References:
42096a0c994SBarry Smith .   1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
42137045ce4SJed Brown 
42237045ce4SJed Brown .seealso: PetscDTLegendreEval()
42337045ce4SJed Brown @*/
42437045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
42537045ce4SJed Brown {
42637045ce4SJed Brown   PetscErrorCode ierr;
42737045ce4SJed Brown   PetscInt       i;
42837045ce4SJed Brown   PetscReal      *work;
42937045ce4SJed Brown   PetscScalar    *Z;
43037045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
43137045ce4SJed Brown 
43237045ce4SJed Brown   PetscFunctionBegin;
4330bfcf5a5SMatthew G. Knepley   ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr);
43437045ce4SJed Brown   /* Set up the Golub-Welsch system */
43537045ce4SJed Brown   for (i=0; i<npoints; i++) {
43637045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
43737045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
43837045ce4SJed Brown   }
439dcca6d9dSJed Brown   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
440c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
44137045ce4SJed Brown   LDZ  = N;
44237045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
4438b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
44437045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
4451c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
44637045ce4SJed Brown 
44737045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
44837045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
44937045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
45019a57d60SBarry Smith     if (x[i] == -0.0) x[i] = 0.0;
45137045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
4520d644c17SKarl Rupp 
45388393a60SJed Brown     w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
45437045ce4SJed Brown   }
45537045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
45637045ce4SJed Brown   PetscFunctionReturn(0);
45737045ce4SJed Brown }
458194825f6SJed Brown 
459744bafbcSMatthew G. Knepley /*@
460744bafbcSMatthew G. Knepley   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
461744bafbcSMatthew G. Knepley 
462744bafbcSMatthew G. Knepley   Not Collective
463744bafbcSMatthew G. Knepley 
464744bafbcSMatthew G. Knepley   Input Arguments:
465744bafbcSMatthew G. Knepley + dim     - The spatial dimension
466744bafbcSMatthew G. Knepley . npoints - number of points in one dimension
467744bafbcSMatthew G. Knepley . a       - left end of interval (often-1)
468744bafbcSMatthew G. Knepley - b       - right end of interval (often +1)
469744bafbcSMatthew G. Knepley 
470744bafbcSMatthew G. Knepley   Output Argument:
471744bafbcSMatthew G. Knepley . q - A PetscQuadrature object
472744bafbcSMatthew G. Knepley 
473744bafbcSMatthew G. Knepley   Level: intermediate
474744bafbcSMatthew G. Knepley 
475744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
476744bafbcSMatthew G. Knepley @*/
477744bafbcSMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
478744bafbcSMatthew G. Knepley {
479744bafbcSMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k;
480744bafbcSMatthew G. Knepley   PetscReal     *x, *w, *xw, *ww;
481744bafbcSMatthew G. Knepley   PetscErrorCode ierr;
482744bafbcSMatthew G. Knepley 
483744bafbcSMatthew G. Knepley   PetscFunctionBegin;
484744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
485744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr);
486744bafbcSMatthew G. Knepley   /* Set up the Golub-Welsch system */
487744bafbcSMatthew G. Knepley   switch (dim) {
488744bafbcSMatthew G. Knepley   case 0:
489744bafbcSMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
490744bafbcSMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
491744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
492744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
493744bafbcSMatthew G. Knepley     x[0] = 0.0;
494744bafbcSMatthew G. Knepley     w[0] = 1.0;
495744bafbcSMatthew G. Knepley     break;
496744bafbcSMatthew G. Knepley   case 1:
497744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr);
498744bafbcSMatthew G. Knepley     break;
499744bafbcSMatthew G. Knepley   case 2:
500744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
501744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
502744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
503744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
504744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+0] = xw[i];
505744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+1] = xw[j];
506744bafbcSMatthew G. Knepley         w[i*npoints+j]         = ww[i] * ww[j];
507744bafbcSMatthew G. Knepley       }
508744bafbcSMatthew G. Knepley     }
509744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
510744bafbcSMatthew G. Knepley     break;
511744bafbcSMatthew G. Knepley   case 3:
512744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
513744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
514744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
515744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
516744bafbcSMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
517744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
518744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
519744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
520744bafbcSMatthew G. Knepley           w[(i*npoints+j)*npoints+k]         = ww[i] * ww[j] * ww[k];
521744bafbcSMatthew G. Knepley         }
522744bafbcSMatthew G. Knepley       }
523744bafbcSMatthew G. Knepley     }
524744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
525744bafbcSMatthew G. Knepley     break;
526744bafbcSMatthew G. Knepley   default:
527744bafbcSMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
528744bafbcSMatthew G. Knepley   }
529744bafbcSMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
530bcede257SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*q, npoints);CHKERRQ(ierr);
531744bafbcSMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr);
532744bafbcSMatthew G. Knepley   PetscFunctionReturn(0);
533744bafbcSMatthew G. Knepley }
534744bafbcSMatthew G. Knepley 
535494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
536494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
537494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
538494e7359SMatthew G. Knepley {
539494e7359SMatthew G. Knepley   PetscReal f = 1.0;
540494e7359SMatthew G. Knepley   PetscInt  i;
541494e7359SMatthew G. Knepley 
542494e7359SMatthew G. Knepley   PetscFunctionBegin;
543494e7359SMatthew G. Knepley   for (i = 1; i < n+1; ++i) f *= i;
544494e7359SMatthew G. Knepley   *factorial = f;
545494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
546494e7359SMatthew G. Knepley }
547494e7359SMatthew G. Knepley 
548494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
549494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
550494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
551494e7359SMatthew G. Knepley {
552494e7359SMatthew G. Knepley   PetscReal apb, pn1, pn2;
553494e7359SMatthew G. Knepley   PetscInt  k;
554494e7359SMatthew G. Knepley 
555494e7359SMatthew G. Knepley   PetscFunctionBegin;
556494e7359SMatthew G. Knepley   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
557494e7359SMatthew G. Knepley   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
558494e7359SMatthew G. Knepley   apb = a + b;
559494e7359SMatthew G. Knepley   pn2 = 1.0;
560494e7359SMatthew G. Knepley   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
561494e7359SMatthew G. Knepley   *P  = 0.0;
562494e7359SMatthew G. Knepley   for (k = 2; k < n+1; ++k) {
563494e7359SMatthew G. Knepley     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
564494e7359SMatthew G. Knepley     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
565494e7359SMatthew G. Knepley     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
566494e7359SMatthew G. Knepley     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
567494e7359SMatthew G. Knepley 
568494e7359SMatthew G. Knepley     a2  = a2 / a1;
569494e7359SMatthew G. Knepley     a3  = a3 / a1;
570494e7359SMatthew G. Knepley     a4  = a4 / a1;
571494e7359SMatthew G. Knepley     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
572494e7359SMatthew G. Knepley     pn2 = pn1;
573494e7359SMatthew G. Knepley     pn1 = *P;
574494e7359SMatthew G. Knepley   }
575494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
576494e7359SMatthew G. Knepley }
577494e7359SMatthew G. Knepley 
578494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
579494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
580494e7359SMatthew G. Knepley {
581494e7359SMatthew G. Knepley   PetscReal      nP;
582494e7359SMatthew G. Knepley   PetscErrorCode ierr;
583494e7359SMatthew G. Knepley 
584494e7359SMatthew G. Knepley   PetscFunctionBegin;
585494e7359SMatthew G. Knepley   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
586494e7359SMatthew G. Knepley   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
587494e7359SMatthew G. Knepley   *P   = 0.5 * (a + b + n + 1) * nP;
588494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
589494e7359SMatthew G. Knepley }
590494e7359SMatthew G. Knepley 
591494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
592494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
593494e7359SMatthew G. Knepley {
594494e7359SMatthew G. Knepley   PetscFunctionBegin;
595494e7359SMatthew G. Knepley   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
596494e7359SMatthew G. Knepley   *eta = y;
597494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
598494e7359SMatthew G. Knepley }
599494e7359SMatthew G. Knepley 
600494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
601494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
602494e7359SMatthew G. Knepley {
603494e7359SMatthew G. Knepley   PetscFunctionBegin;
604494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
605494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
606494e7359SMatthew G. Knepley   *zeta = z;
607494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
608494e7359SMatthew G. Knepley }
609494e7359SMatthew G. Knepley 
610494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
611494e7359SMatthew G. Knepley {
612494e7359SMatthew G. Knepley   PetscInt       maxIter = 100;
613494e7359SMatthew G. Knepley   PetscReal      eps     = 1.0e-8;
614a8291ba1SSatish Balay   PetscReal      a1, a2, a3, a4, a5, a6;
615494e7359SMatthew G. Knepley   PetscInt       k;
616494e7359SMatthew G. Knepley   PetscErrorCode ierr;
617494e7359SMatthew G. Knepley 
618494e7359SMatthew G. Knepley   PetscFunctionBegin;
619a8291ba1SSatish Balay 
6208b49ba18SBarry Smith   a1      = PetscPowReal(2.0, a+b+1);
621a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA)
6220646a658SBarry Smith   a2      = PetscTGamma(a + npoints + 1);
6230646a658SBarry Smith   a3      = PetscTGamma(b + npoints + 1);
6240646a658SBarry Smith   a4      = PetscTGamma(a + b + npoints + 1);
625a8291ba1SSatish Balay #else
62629bcbfd0SToby Isaac   {
627d24bbb91SToby Isaac     PetscInt ia, ib;
62829bcbfd0SToby Isaac 
629d24bbb91SToby Isaac     ia = (PetscInt) a;
630d24bbb91SToby Isaac     ib = (PetscInt) b;
631d24bbb91SToby Isaac     if (ia == a && ib == b && ia + npoints + 1 > 0 && ib + npoints + 1 > 0 && ia + ib + npoints + 1 > 0) { /* All gamma(x) terms are (x-1)! terms */
632d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr);
633d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr);
634d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr);
63529bcbfd0SToby Isaac     } else {
636a8291ba1SSatish Balay       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
63729bcbfd0SToby Isaac     }
63829bcbfd0SToby Isaac   }
639a8291ba1SSatish Balay #endif
640a8291ba1SSatish Balay 
641494e7359SMatthew G. Knepley   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
642494e7359SMatthew G. Knepley   a6   = a1 * a2 * a3 / a4 / a5;
643494e7359SMatthew G. Knepley   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
644494e7359SMatthew G. Knepley    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
645494e7359SMatthew G. Knepley   for (k = 0; k < npoints; ++k) {
6468b49ba18SBarry Smith     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
647494e7359SMatthew G. Knepley     PetscInt  j;
648494e7359SMatthew G. Knepley 
649494e7359SMatthew G. Knepley     if (k > 0) r = 0.5 * (r + x[k-1]);
650494e7359SMatthew G. Knepley     for (j = 0; j < maxIter; ++j) {
651494e7359SMatthew G. Knepley       PetscReal s = 0.0, delta, f, fp;
652494e7359SMatthew G. Knepley       PetscInt  i;
653494e7359SMatthew G. Knepley 
654494e7359SMatthew G. Knepley       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
655494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
656494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
657494e7359SMatthew G. Knepley       delta = f / (fp - f * s);
658494e7359SMatthew G. Knepley       r     = r - delta;
65977b4d14cSPeter Brune       if (PetscAbsReal(delta) < eps) break;
660494e7359SMatthew G. Knepley     }
661494e7359SMatthew G. Knepley     x[k] = r;
662494e7359SMatthew G. Knepley     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
663494e7359SMatthew G. Knepley     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
664494e7359SMatthew G. Knepley   }
665494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
666494e7359SMatthew G. Knepley }
667494e7359SMatthew G. Knepley 
668fd9d31fbSMatthew G. Knepley /*@C
669494e7359SMatthew G. Knepley   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
670494e7359SMatthew G. Knepley 
671494e7359SMatthew G. Knepley   Not Collective
672494e7359SMatthew G. Knepley 
673494e7359SMatthew G. Knepley   Input Arguments:
674494e7359SMatthew G. Knepley + dim   - The simplex dimension
675744bafbcSMatthew G. Knepley . order - The number of points in one dimension
676494e7359SMatthew G. Knepley . a     - left end of interval (often-1)
677494e7359SMatthew G. Knepley - b     - right end of interval (often +1)
678494e7359SMatthew G. Knepley 
679744bafbcSMatthew G. Knepley   Output Argument:
680552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
681494e7359SMatthew G. Knepley 
682494e7359SMatthew G. Knepley   Level: intermediate
683494e7359SMatthew G. Knepley 
684494e7359SMatthew G. Knepley   References:
68596a0c994SBarry Smith .  1. - Karniadakis and Sherwin.  FIAT
686494e7359SMatthew G. Knepley 
687744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
688494e7359SMatthew G. Knepley @*/
689552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
690494e7359SMatthew G. Knepley {
691552aa4f7SMatthew G. Knepley   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
692494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
693494e7359SMatthew G. Knepley   PetscInt       i, j, k;
694494e7359SMatthew G. Knepley   PetscErrorCode ierr;
695494e7359SMatthew G. Knepley 
696494e7359SMatthew G. Knepley   PetscFunctionBegin;
697494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
698785e854fSJed Brown   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
699785e854fSJed Brown   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
700494e7359SMatthew G. Knepley   switch (dim) {
701707aa5c5SMatthew G. Knepley   case 0:
702707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
703707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
704785e854fSJed Brown     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
705785e854fSJed Brown     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
706707aa5c5SMatthew G. Knepley     x[0] = 0.0;
707707aa5c5SMatthew G. Knepley     w[0] = 1.0;
708707aa5c5SMatthew G. Knepley     break;
709494e7359SMatthew G. Knepley   case 1:
710552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
711494e7359SMatthew G. Knepley     break;
712494e7359SMatthew G. Knepley   case 2:
713dcca6d9dSJed Brown     ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr);
714552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
715552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
716552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
717552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
718552aa4f7SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
719552aa4f7SMatthew G. Knepley         w[i*order+j] = 0.5 * wx[i] * wy[j];
720494e7359SMatthew G. Knepley       }
721494e7359SMatthew G. Knepley     }
722494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
723494e7359SMatthew G. Knepley     break;
724494e7359SMatthew G. Knepley   case 3:
725dcca6d9dSJed Brown     ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr);
726552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
727552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
728552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
729552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
730552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
731552aa4f7SMatthew G. Knepley         for (k = 0; k < order; ++k) {
732552aa4f7SMatthew G. Knepley           ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*order+j)*order+k)*3+0], &x[((i*order+j)*order+k)*3+1], &x[((i*order+j)*order+k)*3+2]);CHKERRQ(ierr);
733552aa4f7SMatthew G. Knepley           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
734494e7359SMatthew G. Knepley         }
735494e7359SMatthew G. Knepley       }
736494e7359SMatthew G. Knepley     }
737494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
738494e7359SMatthew G. Knepley     break;
739494e7359SMatthew G. Knepley   default:
740494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
741494e7359SMatthew G. Knepley   }
74221454ff5SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
743bcede257SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*q, order);CHKERRQ(ierr);
74421454ff5SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
745494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
746494e7359SMatthew G. Knepley }
747494e7359SMatthew G. Knepley 
748b3c0f97bSTom Klotz /*@C
749b3c0f97bSTom Klotz   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
750b3c0f97bSTom Klotz 
751b3c0f97bSTom Klotz   Not Collective
752b3c0f97bSTom Klotz 
753b3c0f97bSTom Klotz   Input Arguments:
754b3c0f97bSTom Klotz + dim   - The cell dimension
755b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l
756b3c0f97bSTom Klotz . a     - left end of interval (often-1)
757b3c0f97bSTom Klotz - b     - right end of interval (often +1)
758b3c0f97bSTom Klotz 
759b3c0f97bSTom Klotz   Output Argument:
760b3c0f97bSTom Klotz . q - A PetscQuadrature object
761b3c0f97bSTom Klotz 
762b3c0f97bSTom Klotz   Level: intermediate
763b3c0f97bSTom Klotz 
764b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature()
765b3c0f97bSTom Klotz @*/
766b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
767b3c0f97bSTom Klotz {
768b3c0f97bSTom Klotz   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
769b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;                  /* Half-width of the integration interval */
770b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;                  /* Center of the integration interval */
771b3c0f97bSTom Klotz   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
772d84b4d08SMatthew G. Knepley   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
773b3c0f97bSTom Klotz   PetscReal       wk    = 0.5*PETSC_PI;              /* Quadrature weight at x_k */
774b3c0f97bSTom Klotz   PetscReal      *x, *w;
775b3c0f97bSTom Klotz   PetscInt        K, k, npoints;
776b3c0f97bSTom Klotz   PetscErrorCode  ierr;
777b3c0f97bSTom Klotz 
778b3c0f97bSTom Klotz   PetscFunctionBegin;
779b3c0f97bSTom Klotz   if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
780b3c0f97bSTom Klotz   if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
781b3c0f97bSTom Klotz   /* Find K such that the weights are < 32 digits of precision */
782b3c0f97bSTom Klotz   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
7839add2064SThomas Klotz     wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
784b3c0f97bSTom Klotz   }
785b3c0f97bSTom Klotz   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
786b3c0f97bSTom Klotz   ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr);
787b3c0f97bSTom Klotz   npoints = 2*K-1;
788b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
789b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
790b3c0f97bSTom Klotz   /* Center term */
791b3c0f97bSTom Klotz   x[0] = beta;
792b3c0f97bSTom Klotz   w[0] = 0.5*alpha*PETSC_PI;
793b3c0f97bSTom Klotz   for (k = 1; k < K; ++k) {
7949add2064SThomas Klotz     wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
7959add2064SThomas Klotz     xk = tanh(0.5*PETSC_PI*PetscSinhReal(k*h));
796b3c0f97bSTom Klotz     x[2*k-1] = -alpha*xk+beta;
797b3c0f97bSTom Klotz     w[2*k-1] = wk;
798b3c0f97bSTom Klotz     x[2*k+0] =  alpha*xk+beta;
799b3c0f97bSTom Klotz     w[2*k+0] = wk;
800b3c0f97bSTom Klotz   }
801b3c0f97bSTom Klotz   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
802b3c0f97bSTom Klotz   PetscFunctionReturn(0);
803b3c0f97bSTom Klotz }
804b3c0f97bSTom Klotz 
805b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
806b3c0f97bSTom Klotz {
807b3c0f97bSTom Klotz   const PetscInt  p     = 16;        /* Digits of precision in the evaluation */
808b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;  /* Half-width of the integration interval */
809b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;  /* Center of the integration interval */
810b3c0f97bSTom Klotz   PetscReal       h     = 1.0;       /* Step size, length between x_k */
811b3c0f97bSTom Klotz   PetscInt        l     = 0;         /* Level of refinement, h = 2^{-l} */
812b3c0f97bSTom Klotz   PetscReal       osum  = 0.0;       /* Integral on last level */
813b3c0f97bSTom Klotz   PetscReal       psum  = 0.0;       /* Integral on the level before the last level */
814b3c0f97bSTom Klotz   PetscReal       sum;               /* Integral on current level */
815446c295cSMatthew G. Knepley   PetscReal       yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
816b3c0f97bSTom Klotz   PetscReal       lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
817b3c0f97bSTom Klotz   PetscReal       wk;                /* Quadrature weight at x_k */
818b3c0f97bSTom Klotz   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
819b3c0f97bSTom Klotz   PetscInt        d;                 /* Digits of precision in the integral */
820b3c0f97bSTom Klotz 
821b3c0f97bSTom Klotz   PetscFunctionBegin;
822b3c0f97bSTom Klotz   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
823b3c0f97bSTom Klotz   /* Center term */
824b3c0f97bSTom Klotz   func(beta, &lval);
825b3c0f97bSTom Klotz   sum = 0.5*alpha*PETSC_PI*lval;
826b3c0f97bSTom Klotz   /* */
827b3c0f97bSTom Klotz   do {
828b3c0f97bSTom Klotz     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
829b3c0f97bSTom Klotz     PetscInt  k = 1;
830b3c0f97bSTom Klotz 
831b3c0f97bSTom Klotz     ++l;
832b3c0f97bSTom Klotz     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
833b3c0f97bSTom Klotz     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
834b3c0f97bSTom Klotz     psum = osum;
835b3c0f97bSTom Klotz     osum = sum;
836b3c0f97bSTom Klotz     h   *= 0.5;
837b3c0f97bSTom Klotz     sum *= 0.5;
838b3c0f97bSTom Klotz     do {
8399add2064SThomas Klotz       wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
840446c295cSMatthew G. Knepley       yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
841446c295cSMatthew G. Knepley       lx = -alpha*(1.0 - yk)+beta;
842446c295cSMatthew G. Knepley       rx =  alpha*(1.0 - yk)+beta;
843b3c0f97bSTom Klotz       func(lx, &lval);
844b3c0f97bSTom Klotz       func(rx, &rval);
845b3c0f97bSTom Klotz       lterm   = alpha*wk*lval;
846b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
847b3c0f97bSTom Klotz       sum    += lterm;
848b3c0f97bSTom Klotz       rterm   = alpha*wk*rval;
849b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
850b3c0f97bSTom Klotz       sum    += rterm;
851b3c0f97bSTom Klotz       ++k;
852b3c0f97bSTom Klotz       /* Only need to evaluate every other point on refined levels */
853b3c0f97bSTom Klotz       if (l != 1) ++k;
8549add2064SThomas Klotz     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
855b3c0f97bSTom Klotz 
856b3c0f97bSTom Klotz     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
857b3c0f97bSTom Klotz     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
858b3c0f97bSTom Klotz     d3 = PetscLog10Real(maxTerm) - p;
85909d48545SBarry Smith     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
86009d48545SBarry Smith     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
861b3c0f97bSTom Klotz     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
8629add2064SThomas Klotz   } while (d < digits && l < 12);
863b3c0f97bSTom Klotz   *sol = sum;
864e510cb1fSThomas Klotz 
865b3c0f97bSTom Klotz   PetscFunctionReturn(0);
866b3c0f97bSTom Klotz }
867b3c0f97bSTom Klotz 
86829f144ccSMatthew G. Knepley #ifdef PETSC_HAVE_MPFR
86929f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
87029f144ccSMatthew G. Knepley {
871e510cb1fSThomas Klotz   const PetscInt  safetyFactor = 2;  /* Calculate abcissa until 2*p digits */
87229f144ccSMatthew G. Knepley   PetscInt        l            = 0;  /* Level of refinement, h = 2^{-l} */
87329f144ccSMatthew G. Knepley   mpfr_t          alpha;             /* Half-width of the integration interval */
87429f144ccSMatthew G. Knepley   mpfr_t          beta;              /* Center of the integration interval */
87529f144ccSMatthew G. Knepley   mpfr_t          h;                 /* Step size, length between x_k */
87629f144ccSMatthew G. Knepley   mpfr_t          osum;              /* Integral on last level */
87729f144ccSMatthew G. Knepley   mpfr_t          psum;              /* Integral on the level before the last level */
87829f144ccSMatthew G. Knepley   mpfr_t          sum;               /* Integral on current level */
87929f144ccSMatthew G. Knepley   mpfr_t          yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
88029f144ccSMatthew G. Knepley   mpfr_t          lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
88129f144ccSMatthew G. Knepley   mpfr_t          wk;                /* Quadrature weight at x_k */
88229f144ccSMatthew G. Knepley   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
88329f144ccSMatthew G. Knepley   PetscInt        d;                 /* Digits of precision in the integral */
88429f144ccSMatthew G. Knepley   mpfr_t          pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
88529f144ccSMatthew G. Knepley 
88629f144ccSMatthew G. Knepley   PetscFunctionBegin;
88729f144ccSMatthew G. Knepley   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
88829f144ccSMatthew G. Knepley   /* Create high precision storage */
889c9f744b5SMatthew 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);
89029f144ccSMatthew G. Knepley   /* Initialization */
89129f144ccSMatthew G. Knepley   mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
89229f144ccSMatthew G. Knepley   mpfr_set_d(beta,  0.5*(b+a), MPFR_RNDN);
89329f144ccSMatthew G. Knepley   mpfr_set_d(osum,  0.0,       MPFR_RNDN);
89429f144ccSMatthew G. Knepley   mpfr_set_d(psum,  0.0,       MPFR_RNDN);
89529f144ccSMatthew G. Knepley   mpfr_set_d(h,     1.0,       MPFR_RNDN);
89629f144ccSMatthew G. Knepley   mpfr_const_pi(pi2, MPFR_RNDN);
89729f144ccSMatthew G. Knepley   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
89829f144ccSMatthew G. Knepley   /* Center term */
89929f144ccSMatthew G. Knepley   func(0.5*(b+a), &lval);
90029f144ccSMatthew G. Knepley   mpfr_set(sum, pi2, MPFR_RNDN);
90129f144ccSMatthew G. Knepley   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
90229f144ccSMatthew G. Knepley   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
90329f144ccSMatthew G. Knepley   /* */
90429f144ccSMatthew G. Knepley   do {
90529f144ccSMatthew G. Knepley     PetscReal d1, d2, d3, d4;
90629f144ccSMatthew G. Knepley     PetscInt  k = 1;
90729f144ccSMatthew G. Knepley 
90829f144ccSMatthew G. Knepley     ++l;
90929f144ccSMatthew G. Knepley     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
91029f144ccSMatthew G. Knepley     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
91129f144ccSMatthew G. Knepley     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
91229f144ccSMatthew G. Knepley     mpfr_set(psum, osum, MPFR_RNDN);
91329f144ccSMatthew G. Knepley     mpfr_set(osum,  sum, MPFR_RNDN);
91429f144ccSMatthew G. Knepley     mpfr_mul_d(h,   h,   0.5, MPFR_RNDN);
91529f144ccSMatthew G. Knepley     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
91629f144ccSMatthew G. Knepley     do {
91729f144ccSMatthew G. Knepley       mpfr_set_si(kh, k, MPFR_RNDN);
91829f144ccSMatthew G. Knepley       mpfr_mul(kh, kh, h, MPFR_RNDN);
91929f144ccSMatthew G. Knepley       /* Weight */
92029f144ccSMatthew G. Knepley       mpfr_set(wk, h, MPFR_RNDN);
92129f144ccSMatthew G. Knepley       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
92229f144ccSMatthew G. Knepley       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
92329f144ccSMatthew G. Knepley       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
92429f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
92529f144ccSMatthew G. Knepley       mpfr_sqr(tmp, tmp, MPFR_RNDN);
92629f144ccSMatthew G. Knepley       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
92729f144ccSMatthew G. Knepley       mpfr_div(wk, wk, tmp, MPFR_RNDN);
92829f144ccSMatthew G. Knepley       /* Abscissa */
92929f144ccSMatthew G. Knepley       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
93029f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
93129f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
93229f144ccSMatthew G. Knepley       mpfr_exp(tmp, msinh, MPFR_RNDN);
93329f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
93429f144ccSMatthew G. Knepley       /* Quadrature points */
93529f144ccSMatthew G. Knepley       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
93629f144ccSMatthew G. Knepley       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
93729f144ccSMatthew G. Knepley       mpfr_add(lx, lx, beta, MPFR_RNDU);
93829f144ccSMatthew G. Knepley       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
93929f144ccSMatthew G. Knepley       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
94029f144ccSMatthew G. Knepley       mpfr_add(rx, rx, beta, MPFR_RNDD);
94129f144ccSMatthew G. Knepley       /* Evaluation */
94229f144ccSMatthew G. Knepley       func(mpfr_get_d(lx, MPFR_RNDU), &lval);
94329f144ccSMatthew G. Knepley       func(mpfr_get_d(rx, MPFR_RNDD), &rval);
94429f144ccSMatthew G. Knepley       /* Update */
94529f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
94629f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
94729f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
94829f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
94929f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
95029f144ccSMatthew G. Knepley       mpfr_set(curTerm, tmp, MPFR_RNDN);
95129f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
95229f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
95329f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
95429f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
95529f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
95629f144ccSMatthew G. Knepley       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
95729f144ccSMatthew G. Knepley       ++k;
95829f144ccSMatthew G. Knepley       /* Only need to evaluate every other point on refined levels */
95929f144ccSMatthew G. Knepley       if (l != 1) ++k;
96029f144ccSMatthew G. Knepley       mpfr_log10(tmp, wk, MPFR_RNDN);
96129f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
962c9f744b5SMatthew G. Knepley     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
96329f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
96429f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
96529f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
96629f144ccSMatthew G. Knepley     d1 = mpfr_get_d(tmp, MPFR_RNDN);
96729f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
96829f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
96929f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
97029f144ccSMatthew G. Knepley     d2 = mpfr_get_d(tmp, MPFR_RNDN);
97129f144ccSMatthew G. Knepley     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
972c9f744b5SMatthew G. Knepley     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
97329f144ccSMatthew G. Knepley     mpfr_log10(tmp, curTerm, MPFR_RNDN);
97429f144ccSMatthew G. Knepley     d4 = mpfr_get_d(tmp, MPFR_RNDN);
97529f144ccSMatthew G. Knepley     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
976b0649871SThomas Klotz   } while (d < digits && l < 8);
97729f144ccSMatthew G. Knepley   *sol = mpfr_get_d(sum, MPFR_RNDN);
97829f144ccSMatthew G. Knepley   /* Cleanup */
97929f144ccSMatthew G. Knepley   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
98029f144ccSMatthew G. Knepley   PetscFunctionReturn(0);
98129f144ccSMatthew G. Knepley }
982d525116cSMatthew G. Knepley #else
983*fbfcfee5SBarry Smith 
984d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
985d525116cSMatthew G. Knepley {
986d525116cSMatthew G. Knepley   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
987d525116cSMatthew G. Knepley }
98829f144ccSMatthew G. Knepley #endif
98929f144ccSMatthew G. Knepley 
990194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
991194825f6SJed Brown  * A in column-major format
992194825f6SJed Brown  * Ainv in row-major format
993194825f6SJed Brown  * tau has length m
994194825f6SJed Brown  * worksize must be >= max(1,n)
995194825f6SJed Brown  */
996194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
997194825f6SJed Brown {
998194825f6SJed Brown   PetscErrorCode ierr;
999194825f6SJed Brown   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
1000194825f6SJed Brown   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
1001194825f6SJed Brown 
1002194825f6SJed Brown   PetscFunctionBegin;
1003194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
1004194825f6SJed Brown   {
1005194825f6SJed Brown     PetscInt i,j;
1006dcca6d9dSJed Brown     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
1007194825f6SJed Brown     for (j=0; j<n; j++) {
1008194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1009194825f6SJed Brown     }
1010194825f6SJed Brown     mstride = m;
1011194825f6SJed Brown   }
1012194825f6SJed Brown #else
1013194825f6SJed Brown   A = A_in;
1014194825f6SJed Brown   Ainv = Ainv_out;
1015194825f6SJed Brown #endif
1016194825f6SJed Brown 
1017194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
1018194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
1019194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
1020194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
1021194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1022001a771dSBarry Smith   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1023194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1024194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1025194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1026194825f6SJed Brown 
1027194825f6SJed Brown   /* Extract an explicit representation of Q */
1028194825f6SJed Brown   Q = Ainv;
1029194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
1030194825f6SJed Brown   K = N;                        /* full rank */
1031001a771dSBarry Smith   PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1032194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1033194825f6SJed Brown 
1034194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1035194825f6SJed Brown   Alpha = 1.0;
1036194825f6SJed Brown   ldb = lda;
1037001a771dSBarry Smith   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1038194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
1039194825f6SJed Brown 
1040194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
1041194825f6SJed Brown   {
1042194825f6SJed Brown     PetscInt i;
1043194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1044194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
1045194825f6SJed Brown   }
1046194825f6SJed Brown #endif
1047194825f6SJed Brown   PetscFunctionReturn(0);
1048194825f6SJed Brown }
1049194825f6SJed Brown 
1050194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1051194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1052194825f6SJed Brown {
1053194825f6SJed Brown   PetscErrorCode ierr;
1054194825f6SJed Brown   PetscReal      *Bv;
1055194825f6SJed Brown   PetscInt       i,j;
1056194825f6SJed Brown 
1057194825f6SJed Brown   PetscFunctionBegin;
1058785e854fSJed Brown   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
1059194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
1060194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
1061194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1062194825f6SJed Brown   for (i=0; i<ninterval; i++) {
1063194825f6SJed Brown     for (j=0; j<ndegree; j++) {
1064194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1065194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1066194825f6SJed Brown     }
1067194825f6SJed Brown   }
1068194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
1069194825f6SJed Brown   PetscFunctionReturn(0);
1070194825f6SJed Brown }
1071194825f6SJed Brown 
1072194825f6SJed Brown /*@
1073194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
1074194825f6SJed Brown 
1075194825f6SJed Brown    Not Collective
1076194825f6SJed Brown 
1077194825f6SJed Brown    Input Arguments:
1078194825f6SJed Brown +  degree - degree of reconstruction polynomial
1079194825f6SJed Brown .  nsource - number of source intervals
1080194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1081194825f6SJed Brown .  ntarget - number of target intervals
1082194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
1083194825f6SJed Brown 
1084194825f6SJed Brown    Output Arguments:
1085194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
1086194825f6SJed Brown 
1087194825f6SJed Brown    Level: advanced
1088194825f6SJed Brown 
1089194825f6SJed Brown .seealso: PetscDTLegendreEval()
1090194825f6SJed Brown @*/
1091194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1092194825f6SJed Brown {
1093194825f6SJed Brown   PetscErrorCode ierr;
1094194825f6SJed Brown   PetscInt       i,j,k,*bdegrees,worksize;
1095194825f6SJed Brown   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1096194825f6SJed Brown   PetscScalar    *tau,*work;
1097194825f6SJed Brown 
1098194825f6SJed Brown   PetscFunctionBegin;
1099194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
1100194825f6SJed Brown   PetscValidRealPointer(targetx,5);
1101194825f6SJed Brown   PetscValidRealPointer(R,6);
1102194825f6SJed Brown   if (degree >= nsource) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Reconstruction degree %D must be less than number of source intervals %D",degree,nsource);
1103194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
1104194825f6SJed Brown   for (i=0; i<nsource; i++) {
110557622a8eSBarry Smith     if (sourcex[i] >= sourcex[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Source interval %D has negative orientation (%g,%g)",i,(double)sourcex[i],(double)sourcex[i+1]);
1106194825f6SJed Brown   }
1107194825f6SJed Brown   for (i=0; i<ntarget; i++) {
110857622a8eSBarry Smith     if (targetx[i] >= targetx[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Target interval %D has negative orientation (%g,%g)",i,(double)targetx[i],(double)targetx[i+1]);
1109194825f6SJed Brown   }
1110194825f6SJed Brown #endif
1111194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
1112194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1113194825f6SJed Brown   center = (xmin + xmax)/2;
1114194825f6SJed Brown   hscale = (xmax - xmin)/2;
1115194825f6SJed Brown   worksize = nsource;
1116dcca6d9dSJed Brown   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
1117dcca6d9dSJed Brown   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
1118194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1119194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1120194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
1121194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
1122194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1123194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
1124194825f6SJed Brown   for (i=0; i<ntarget; i++) {
1125194825f6SJed Brown     PetscReal rowsum = 0;
1126194825f6SJed Brown     for (j=0; j<nsource; j++) {
1127194825f6SJed Brown       PetscReal sum = 0;
1128194825f6SJed Brown       for (k=0; k<degree+1; k++) {
1129194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1130194825f6SJed Brown       }
1131194825f6SJed Brown       R[i*nsource+j] = sum;
1132194825f6SJed Brown       rowsum += sum;
1133194825f6SJed Brown     }
1134194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1135194825f6SJed Brown   }
1136194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
1137194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
1138194825f6SJed Brown   PetscFunctionReturn(0);
1139194825f6SJed Brown }
1140