xref: /petsc/src/dm/dt/interface/dt.c (revision 9add2064e55d1bd249b32f71b810cdcea18c3ce2)
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
7a6fc04d9SSatish Balay 
80c35b76eSJed Brown #include <petscdt.h>            /*I "petscdt.h" I*/
937045ce4SJed Brown #include <petscblaslapack.h>
10af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
11af0996ceSBarry Smith #include <petsc/private/dtimpl.h>
12665c2dedSJed Brown #include <petscviewer.h>
1359804f93SMatthew G. Knepley #include <petscdmplex.h>
1459804f93SMatthew G. Knepley #include <petscdmshell.h>
1537045ce4SJed Brown 
160bfcf5a5SMatthew G. Knepley static PetscBool GaussCite       = PETSC_FALSE;
170bfcf5a5SMatthew G. Knepley const char       GaussCitation[] = "@article{GolubWelsch1969,\n"
180bfcf5a5SMatthew G. Knepley                                    "  author  = {Golub and Welsch},\n"
190bfcf5a5SMatthew G. Knepley                                    "  title   = {Calculation of Quadrature Rules},\n"
200bfcf5a5SMatthew G. Knepley                                    "  journal = {Math. Comp.},\n"
210bfcf5a5SMatthew G. Knepley                                    "  volume  = {23},\n"
220bfcf5a5SMatthew G. Knepley                                    "  number  = {106},\n"
230bfcf5a5SMatthew G. Knepley                                    "  pages   = {221--230},\n"
240bfcf5a5SMatthew G. Knepley                                    "  year    = {1969}\n}\n";
250bfcf5a5SMatthew G. Knepley 
2637045ce4SJed Brown #undef __FUNCT__
2721454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureCreate"
2840d8ff71SMatthew G. Knepley /*@
2940d8ff71SMatthew G. Knepley   PetscQuadratureCreate - Create a PetscQuadrature object
3040d8ff71SMatthew G. Knepley 
3140d8ff71SMatthew G. Knepley   Collective on MPI_Comm
3240d8ff71SMatthew G. Knepley 
3340d8ff71SMatthew G. Knepley   Input Parameter:
3440d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object
3540d8ff71SMatthew G. Knepley 
3640d8ff71SMatthew G. Knepley   Output Parameter:
3740d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
3840d8ff71SMatthew G. Knepley 
3940d8ff71SMatthew G. Knepley   Level: beginner
4040d8ff71SMatthew G. Knepley 
4140d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, create
4240d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
4340d8ff71SMatthew G. Knepley @*/
4421454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
4521454ff5SMatthew G. Knepley {
4621454ff5SMatthew G. Knepley   PetscErrorCode ierr;
4721454ff5SMatthew G. Knepley 
4821454ff5SMatthew G. Knepley   PetscFunctionBegin;
4921454ff5SMatthew G. Knepley   PetscValidPointer(q, 2);
5021454ff5SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
5173107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*q,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
5221454ff5SMatthew G. Knepley   (*q)->dim       = -1;
53bcede257SMatthew G. Knepley   (*q)->order     = -1;
5421454ff5SMatthew G. Knepley   (*q)->numPoints = 0;
5521454ff5SMatthew G. Knepley   (*q)->points    = NULL;
5621454ff5SMatthew G. Knepley   (*q)->weights   = NULL;
5721454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
5821454ff5SMatthew G. Knepley }
5921454ff5SMatthew G. Knepley 
6021454ff5SMatthew G. Knepley #undef __FUNCT__
61c9638911SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDuplicate"
62c9638911SMatthew G. Knepley /*@
63c9638911SMatthew G. Knepley   PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object
64c9638911SMatthew G. Knepley 
65c9638911SMatthew G. Knepley   Collective on PetscQuadrature
66c9638911SMatthew G. Knepley 
67c9638911SMatthew G. Knepley   Input Parameter:
68c9638911SMatthew G. Knepley . q  - The PetscQuadrature object
69c9638911SMatthew G. Knepley 
70c9638911SMatthew G. Knepley   Output Parameter:
71c9638911SMatthew G. Knepley . r  - The new PetscQuadrature object
72c9638911SMatthew G. Knepley 
73c9638911SMatthew G. Knepley   Level: beginner
74c9638911SMatthew G. Knepley 
75c9638911SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, clone
76c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData()
77c9638911SMatthew G. Knepley @*/
78c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
79c9638911SMatthew G. Knepley {
80c9638911SMatthew G. Knepley   PetscInt         order, dim, Nq;
81c9638911SMatthew G. Knepley   const PetscReal *points, *weights;
82c9638911SMatthew G. Knepley   PetscReal       *p, *w;
83c9638911SMatthew G. Knepley   PetscErrorCode   ierr;
84c9638911SMatthew G. Knepley 
85c9638911SMatthew G. Knepley   PetscFunctionBegin;
86c9638911SMatthew G. Knepley   PetscValidPointer(q, 2);
87c9638911SMatthew G. Knepley   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr);
88c9638911SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
89c9638911SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr);
90c9638911SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &Nq, &points, &weights);CHKERRQ(ierr);
91c9638911SMatthew G. Knepley   ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr);
92c9638911SMatthew G. Knepley   ierr = PetscMalloc1(Nq, &w);CHKERRQ(ierr);
93c9638911SMatthew G. Knepley   ierr = PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));CHKERRQ(ierr);
94c9638911SMatthew G. Knepley   ierr = PetscMemcpy(w, weights, Nq * sizeof(PetscReal));CHKERRQ(ierr);
95c9638911SMatthew G. Knepley   ierr = PetscQuadratureSetData(*r, dim, Nq, p, w);CHKERRQ(ierr);
96c9638911SMatthew G. Knepley   PetscFunctionReturn(0);
97c9638911SMatthew G. Knepley }
98c9638911SMatthew G. Knepley 
99c9638911SMatthew G. Knepley #undef __FUNCT__
100bfa639d9SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDestroy"
10140d8ff71SMatthew G. Knepley /*@
10240d8ff71SMatthew G. Knepley   PetscQuadratureDestroy - Destroys a PetscQuadrature object
10340d8ff71SMatthew G. Knepley 
10440d8ff71SMatthew G. Knepley   Collective on PetscQuadrature
10540d8ff71SMatthew G. Knepley 
10640d8ff71SMatthew G. Knepley   Input Parameter:
10740d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
10840d8ff71SMatthew G. Knepley 
10940d8ff71SMatthew G. Knepley   Level: beginner
11040d8ff71SMatthew G. Knepley 
11140d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, destroy
11240d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
11340d8ff71SMatthew G. Knepley @*/
114bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
115bfa639d9SMatthew G. Knepley {
116bfa639d9SMatthew G. Knepley   PetscErrorCode ierr;
117bfa639d9SMatthew G. Knepley 
118bfa639d9SMatthew G. Knepley   PetscFunctionBegin;
11921454ff5SMatthew G. Knepley   if (!*q) PetscFunctionReturn(0);
12021454ff5SMatthew G. Knepley   PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1);
12121454ff5SMatthew G. Knepley   if (--((PetscObject)(*q))->refct > 0) {
12221454ff5SMatthew G. Knepley     *q = NULL;
12321454ff5SMatthew G. Knepley     PetscFunctionReturn(0);
12421454ff5SMatthew G. Knepley   }
12521454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
12621454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
12721454ff5SMatthew G. Knepley   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
12821454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
12921454ff5SMatthew G. Knepley }
13021454ff5SMatthew G. Knepley 
13121454ff5SMatthew G. Knepley #undef __FUNCT__
132bcede257SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetOrder"
133bcede257SMatthew G. Knepley /*@
134bcede257SMatthew G. Knepley   PetscQuadratureGetOrder - Return the quadrature information
135bcede257SMatthew G. Knepley 
136bcede257SMatthew G. Knepley   Not collective
137bcede257SMatthew G. Knepley 
138bcede257SMatthew G. Knepley   Input Parameter:
139bcede257SMatthew G. Knepley . q - The PetscQuadrature object
140bcede257SMatthew G. Knepley 
141bcede257SMatthew G. Knepley   Output Parameter:
142bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
143bcede257SMatthew G. Knepley 
144bcede257SMatthew G. Knepley   Output Parameter:
145bcede257SMatthew G. Knepley 
146bcede257SMatthew G. Knepley   Level: intermediate
147bcede257SMatthew G. Knepley 
148bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
149bcede257SMatthew G. Knepley @*/
150bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
151bcede257SMatthew G. Knepley {
152bcede257SMatthew G. Knepley   PetscFunctionBegin;
153bcede257SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
154bcede257SMatthew G. Knepley   PetscValidPointer(order, 2);
155bcede257SMatthew G. Knepley   *order = q->order;
156bcede257SMatthew G. Knepley   PetscFunctionReturn(0);
157bcede257SMatthew G. Knepley }
158bcede257SMatthew G. Knepley 
159bcede257SMatthew G. Knepley #undef __FUNCT__
160bcede257SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetOrder"
161bcede257SMatthew G. Knepley /*@
162bcede257SMatthew G. Knepley   PetscQuadratureSetOrder - Return the quadrature information
163bcede257SMatthew G. Knepley 
164bcede257SMatthew G. Knepley   Not collective
165bcede257SMatthew G. Knepley 
166bcede257SMatthew G. Knepley   Input Parameters:
167bcede257SMatthew G. Knepley + 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 
172bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
173bcede257SMatthew G. Knepley @*/
174bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
175bcede257SMatthew G. Knepley {
176bcede257SMatthew G. Knepley   PetscFunctionBegin;
177bcede257SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
178bcede257SMatthew G. Knepley   q->order = order;
179bcede257SMatthew G. Knepley   PetscFunctionReturn(0);
180bcede257SMatthew G. Knepley }
181bcede257SMatthew G. Knepley 
182bcede257SMatthew G. Knepley #undef __FUNCT__
18321454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetData"
18440d8ff71SMatthew G. Knepley /*@C
18540d8ff71SMatthew G. Knepley   PetscQuadratureGetData - Returns the data defining the quadrature
18640d8ff71SMatthew G. Knepley 
18740d8ff71SMatthew G. Knepley   Not collective
18840d8ff71SMatthew G. Knepley 
18940d8ff71SMatthew G. Knepley   Input Parameter:
19040d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
19140d8ff71SMatthew G. Knepley 
19240d8ff71SMatthew G. Knepley   Output Parameters:
19340d8ff71SMatthew G. Knepley + dim - The spatial dimension
19440d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
19540d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
19640d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
19740d8ff71SMatthew G. Knepley 
19840d8ff71SMatthew G. Knepley   Level: intermediate
19940d8ff71SMatthew G. Knepley 
20040d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature
20140d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
20240d8ff71SMatthew G. Knepley @*/
20321454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
20421454ff5SMatthew G. Knepley {
20521454ff5SMatthew G. Knepley   PetscFunctionBegin;
20621454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
20721454ff5SMatthew G. Knepley   if (dim) {
20821454ff5SMatthew G. Knepley     PetscValidPointer(dim, 2);
20921454ff5SMatthew G. Knepley     *dim = q->dim;
21021454ff5SMatthew G. Knepley   }
21121454ff5SMatthew G. Knepley   if (npoints) {
21221454ff5SMatthew G. Knepley     PetscValidPointer(npoints, 3);
21321454ff5SMatthew G. Knepley     *npoints = q->numPoints;
21421454ff5SMatthew G. Knepley   }
21521454ff5SMatthew G. Knepley   if (points) {
21621454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
21721454ff5SMatthew G. Knepley     *points = q->points;
21821454ff5SMatthew G. Knepley   }
21921454ff5SMatthew G. Knepley   if (weights) {
22021454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
22121454ff5SMatthew G. Knepley     *weights = q->weights;
22221454ff5SMatthew G. Knepley   }
22321454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
22421454ff5SMatthew G. Knepley }
22521454ff5SMatthew G. Knepley 
22621454ff5SMatthew G. Knepley #undef __FUNCT__
22721454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetData"
22840d8ff71SMatthew G. Knepley /*@C
22940d8ff71SMatthew G. Knepley   PetscQuadratureSetData - Sets the data defining the quadrature
23040d8ff71SMatthew G. Knepley 
23140d8ff71SMatthew G. Knepley   Not collective
23240d8ff71SMatthew G. Knepley 
23340d8ff71SMatthew G. Knepley   Input Parameters:
23440d8ff71SMatthew G. Knepley + q  - The PetscQuadrature object
23540d8ff71SMatthew G. Knepley . dim - The spatial dimension
23640d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
23740d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
23840d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
23940d8ff71SMatthew G. Knepley 
24040d8ff71SMatthew G. Knepley   Level: intermediate
24140d8ff71SMatthew G. Knepley 
24240d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature
24340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
24440d8ff71SMatthew G. Knepley @*/
24521454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
24621454ff5SMatthew G. Knepley {
24721454ff5SMatthew G. Knepley   PetscFunctionBegin;
24821454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
24921454ff5SMatthew G. Knepley   if (dim >= 0)     q->dim       = dim;
25021454ff5SMatthew G. Knepley   if (npoints >= 0) q->numPoints = npoints;
25121454ff5SMatthew G. Knepley   if (points) {
25221454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
25321454ff5SMatthew G. Knepley     q->points = points;
25421454ff5SMatthew G. Knepley   }
25521454ff5SMatthew G. Knepley   if (weights) {
25621454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
25721454ff5SMatthew G. Knepley     q->weights = weights;
25821454ff5SMatthew G. Knepley   }
259f9fd7fdbSMatthew G. Knepley   PetscFunctionReturn(0);
260f9fd7fdbSMatthew G. Knepley }
261f9fd7fdbSMatthew G. Knepley 
262f9fd7fdbSMatthew G. Knepley #undef __FUNCT__
263f9fd7fdbSMatthew G. Knepley #define __FUNCT__ "PetscQuadratureView"
26440d8ff71SMatthew G. Knepley /*@C
26540d8ff71SMatthew G. Knepley   PetscQuadratureView - Views a PetscQuadrature object
26640d8ff71SMatthew G. Knepley 
26740d8ff71SMatthew G. Knepley   Collective on PetscQuadrature
26840d8ff71SMatthew G. Knepley 
26940d8ff71SMatthew G. Knepley   Input Parameters:
27040d8ff71SMatthew G. Knepley + q  - The PetscQuadrature object
27140d8ff71SMatthew G. Knepley - viewer - The PetscViewer object
27240d8ff71SMatthew G. Knepley 
27340d8ff71SMatthew G. Knepley   Level: beginner
27440d8ff71SMatthew G. Knepley 
27540d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, view
27640d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
27740d8ff71SMatthew G. Knepley @*/
278f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
279f9fd7fdbSMatthew G. Knepley {
280f9fd7fdbSMatthew G. Knepley   PetscInt       q, d;
281f9fd7fdbSMatthew G. Knepley   PetscErrorCode ierr;
282f9fd7fdbSMatthew G. Knepley 
283f9fd7fdbSMatthew G. Knepley   PetscFunctionBegin;
28498c3331eSBarry Smith   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr);
28521454ff5SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n  (", quad->numPoints);CHKERRQ(ierr);
28621454ff5SMatthew G. Knepley   for (q = 0; q < quad->numPoints; ++q) {
28721454ff5SMatthew G. Knepley     for (d = 0; d < quad->dim; ++d) {
288f9fd7fdbSMatthew G. Knepley       if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);
289ab15ae43SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
290f9fd7fdbSMatthew G. Knepley     }
291ab15ae43SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr);
292f9fd7fdbSMatthew G. Knepley   }
293bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
294bfa639d9SMatthew G. Knepley }
295bfa639d9SMatthew G. Knepley 
296bfa639d9SMatthew G. Knepley #undef __FUNCT__
29789710940SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureExpandComposite"
29889710940SMatthew G. Knepley /*@C
29989710940SMatthew G. Knepley   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
30089710940SMatthew G. Knepley 
30189710940SMatthew G. Knepley   Not collective
30289710940SMatthew G. Knepley 
30389710940SMatthew G. Knepley   Input Parameter:
30489710940SMatthew G. Knepley + q - The original PetscQuadrature
30589710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into
30689710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement
30789710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement
30889710940SMatthew G. Knepley 
30989710940SMatthew G. Knepley   Output Parameters:
31089710940SMatthew G. Knepley . dim - The dimension
31189710940SMatthew G. Knepley 
31289710940SMatthew G. Knepley   Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
31389710940SMatthew G. Knepley 
31489710940SMatthew G. Knepley   Level: intermediate
31589710940SMatthew G. Knepley 
31689710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
31789710940SMatthew G. Knepley @*/
31889710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
31989710940SMatthew G. Knepley {
32089710940SMatthew G. Knepley   const PetscReal *points,    *weights;
32189710940SMatthew G. Knepley   PetscReal       *pointsRef, *weightsRef;
32289710940SMatthew G. Knepley   PetscInt         dim, order, npoints, npointsRef, c, p, d, e;
32389710940SMatthew G. Knepley   PetscErrorCode   ierr;
32489710940SMatthew G. Knepley 
32589710940SMatthew G. Knepley   PetscFunctionBegin;
32689710940SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
32789710940SMatthew G. Knepley   PetscValidPointer(v0, 3);
32889710940SMatthew G. Knepley   PetscValidPointer(jac, 4);
32989710940SMatthew G. Knepley   PetscValidPointer(qref, 5);
33089710940SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr);
33189710940SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
33289710940SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &npoints, &points, &weights);CHKERRQ(ierr);
33389710940SMatthew G. Knepley   npointsRef = npoints*numSubelements;
33489710940SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr);
33589710940SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef,&weightsRef);CHKERRQ(ierr);
33689710940SMatthew G. Knepley   for (c = 0; c < numSubelements; ++c) {
33789710940SMatthew G. Knepley     for (p = 0; p < npoints; ++p) {
33889710940SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
33989710940SMatthew G. Knepley         pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
34089710940SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
34189710940SMatthew G. Knepley           pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
34289710940SMatthew G. Knepley         }
34389710940SMatthew G. Knepley       }
34489710940SMatthew G. Knepley       /* Could also use detJ here */
34589710940SMatthew G. Knepley       weightsRef[c*npoints+p] = weights[p]/numSubelements;
34689710940SMatthew G. Knepley     }
34789710940SMatthew G. Knepley   }
34889710940SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr);
34989710940SMatthew G. Knepley   ierr = PetscQuadratureSetData(*qref, dim, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr);
35089710940SMatthew G. Knepley   PetscFunctionReturn(0);
35189710940SMatthew G. Knepley }
35289710940SMatthew G. Knepley 
35389710940SMatthew G. Knepley #undef __FUNCT__
35437045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval"
35537045ce4SJed Brown /*@
35637045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
35737045ce4SJed Brown 
35837045ce4SJed Brown    Not Collective
35937045ce4SJed Brown 
36037045ce4SJed Brown    Input Arguments:
36137045ce4SJed Brown +  npoints - number of spatial points to evaluate at
36237045ce4SJed Brown .  points - array of locations to evaluate at
36337045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
36437045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
36537045ce4SJed Brown 
36637045ce4SJed Brown    Output Arguments:
3670298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
3680298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
3690298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
37037045ce4SJed Brown 
37137045ce4SJed Brown    Level: intermediate
37237045ce4SJed Brown 
37337045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
37437045ce4SJed Brown @*/
37537045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
37637045ce4SJed Brown {
37737045ce4SJed Brown   PetscInt i,maxdegree;
37837045ce4SJed Brown 
37937045ce4SJed Brown   PetscFunctionBegin;
38037045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
38137045ce4SJed Brown   maxdegree = degrees[ndegree-1];
38237045ce4SJed Brown   for (i=0; i<npoints; i++) {
38337045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
38437045ce4SJed Brown     PetscInt  j,k;
38537045ce4SJed Brown     x    = points[i];
38637045ce4SJed Brown     pm2  = 0;
38737045ce4SJed Brown     pm1  = 1;
38837045ce4SJed Brown     pd2  = 0;
38937045ce4SJed Brown     pd1  = 0;
39037045ce4SJed Brown     pdd2 = 0;
39137045ce4SJed Brown     pdd1 = 0;
39237045ce4SJed Brown     k    = 0;
39337045ce4SJed Brown     if (degrees[k] == 0) {
39437045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
39537045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
39637045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
39737045ce4SJed Brown       k++;
39837045ce4SJed Brown     }
39937045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
40037045ce4SJed Brown       PetscReal p,d,dd;
40137045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
40237045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
40337045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
40437045ce4SJed Brown       pm2  = pm1;
40537045ce4SJed Brown       pm1  = p;
40637045ce4SJed Brown       pd2  = pd1;
40737045ce4SJed Brown       pd1  = d;
40837045ce4SJed Brown       pdd2 = pdd1;
40937045ce4SJed Brown       pdd1 = dd;
41037045ce4SJed Brown       if (degrees[k] == j) {
41137045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
41237045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
41337045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
41437045ce4SJed Brown       }
41537045ce4SJed Brown     }
41637045ce4SJed Brown   }
41737045ce4SJed Brown   PetscFunctionReturn(0);
41837045ce4SJed Brown }
41937045ce4SJed Brown 
42037045ce4SJed Brown #undef __FUNCT__
42137045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature"
42237045ce4SJed Brown /*@
42337045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
42437045ce4SJed Brown 
42537045ce4SJed Brown    Not Collective
42637045ce4SJed Brown 
42737045ce4SJed Brown    Input Arguments:
42837045ce4SJed Brown +  npoints - number of points
42937045ce4SJed Brown .  a - left end of interval (often-1)
43037045ce4SJed Brown -  b - right end of interval (often +1)
43137045ce4SJed Brown 
43237045ce4SJed Brown    Output Arguments:
43337045ce4SJed Brown +  x - quadrature points
43437045ce4SJed Brown -  w - quadrature weights
43537045ce4SJed Brown 
43637045ce4SJed Brown    Level: intermediate
43737045ce4SJed Brown 
43837045ce4SJed Brown    References:
43937045ce4SJed Brown    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
44037045ce4SJed Brown 
44137045ce4SJed Brown .seealso: PetscDTLegendreEval()
44237045ce4SJed Brown @*/
44337045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
44437045ce4SJed Brown {
44537045ce4SJed Brown   PetscErrorCode ierr;
44637045ce4SJed Brown   PetscInt       i;
44737045ce4SJed Brown   PetscReal      *work;
44837045ce4SJed Brown   PetscScalar    *Z;
44937045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
45037045ce4SJed Brown 
45137045ce4SJed Brown   PetscFunctionBegin;
4520bfcf5a5SMatthew G. Knepley   ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr);
45337045ce4SJed Brown   /* Set up the Golub-Welsch system */
45437045ce4SJed Brown   for (i=0; i<npoints; i++) {
45537045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
45637045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
45737045ce4SJed Brown   }
458dcca6d9dSJed Brown   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
459c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
46037045ce4SJed Brown   LDZ  = N;
46137045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
4628b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
46337045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
4641c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
46537045ce4SJed Brown 
46637045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
46737045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
46837045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
46937045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
4700d644c17SKarl Rupp 
47188393a60SJed Brown     w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
47237045ce4SJed Brown   }
47337045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
47437045ce4SJed Brown   PetscFunctionReturn(0);
47537045ce4SJed Brown }
476194825f6SJed Brown 
477194825f6SJed Brown #undef __FUNCT__
478744bafbcSMatthew G. Knepley #define __FUNCT__ "PetscDTGaussTensorQuadrature"
479744bafbcSMatthew G. Knepley /*@
480744bafbcSMatthew G. Knepley   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
481744bafbcSMatthew G. Knepley 
482744bafbcSMatthew G. Knepley   Not Collective
483744bafbcSMatthew G. Knepley 
484744bafbcSMatthew G. Knepley   Input Arguments:
485744bafbcSMatthew G. Knepley + dim     - The spatial dimension
486744bafbcSMatthew G. Knepley . npoints - number of points in one dimension
487744bafbcSMatthew G. Knepley . a       - left end of interval (often-1)
488744bafbcSMatthew G. Knepley - b       - right end of interval (often +1)
489744bafbcSMatthew G. Knepley 
490744bafbcSMatthew G. Knepley   Output Argument:
491744bafbcSMatthew G. Knepley . q - A PetscQuadrature object
492744bafbcSMatthew G. Knepley 
493744bafbcSMatthew G. Knepley   Level: intermediate
494744bafbcSMatthew G. Knepley 
495744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
496744bafbcSMatthew G. Knepley @*/
497744bafbcSMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
498744bafbcSMatthew G. Knepley {
499744bafbcSMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k;
500744bafbcSMatthew G. Knepley   PetscReal     *x, *w, *xw, *ww;
501744bafbcSMatthew G. Knepley   PetscErrorCode ierr;
502744bafbcSMatthew G. Knepley 
503744bafbcSMatthew G. Knepley   PetscFunctionBegin;
504744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
505744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr);
506744bafbcSMatthew G. Knepley   /* Set up the Golub-Welsch system */
507744bafbcSMatthew G. Knepley   switch (dim) {
508744bafbcSMatthew G. Knepley   case 0:
509744bafbcSMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
510744bafbcSMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
511744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
512744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
513744bafbcSMatthew G. Knepley     x[0] = 0.0;
514744bafbcSMatthew G. Knepley     w[0] = 1.0;
515744bafbcSMatthew G. Knepley     break;
516744bafbcSMatthew G. Knepley   case 1:
517744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr);
518744bafbcSMatthew G. Knepley     break;
519744bafbcSMatthew G. Knepley   case 2:
520744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
521744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
522744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
523744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
524744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+0] = xw[i];
525744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+1] = xw[j];
526744bafbcSMatthew G. Knepley         w[i*npoints+j]         = ww[i] * ww[j];
527744bafbcSMatthew G. Knepley       }
528744bafbcSMatthew G. Knepley     }
529744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
530744bafbcSMatthew G. Knepley     break;
531744bafbcSMatthew G. Knepley   case 3:
532744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
533744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
534744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
535744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
536744bafbcSMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
537744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
538744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
539744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
540744bafbcSMatthew G. Knepley           w[(i*npoints+j)*npoints+k]         = ww[i] * ww[j] * ww[k];
541744bafbcSMatthew G. Knepley         }
542744bafbcSMatthew G. Knepley       }
543744bafbcSMatthew G. Knepley     }
544744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
545744bafbcSMatthew G. Knepley     break;
546744bafbcSMatthew G. Knepley   default:
547744bafbcSMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
548744bafbcSMatthew G. Knepley   }
549744bafbcSMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
550bcede257SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*q, npoints);CHKERRQ(ierr);
551744bafbcSMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr);
552744bafbcSMatthew G. Knepley   PetscFunctionReturn(0);
553744bafbcSMatthew G. Knepley }
554744bafbcSMatthew G. Knepley 
555744bafbcSMatthew G. Knepley #undef __FUNCT__
556494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal"
557494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
558494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
559494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
560494e7359SMatthew G. Knepley {
561494e7359SMatthew G. Knepley   PetscReal f = 1.0;
562494e7359SMatthew G. Knepley   PetscInt  i;
563494e7359SMatthew G. Knepley 
564494e7359SMatthew G. Knepley   PetscFunctionBegin;
565494e7359SMatthew G. Knepley   for (i = 1; i < n+1; ++i) f *= i;
566494e7359SMatthew G. Knepley   *factorial = f;
567494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
568494e7359SMatthew G. Knepley }
569494e7359SMatthew G. Knepley 
570494e7359SMatthew G. Knepley #undef __FUNCT__
571494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi"
572494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
573494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
574494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
575494e7359SMatthew G. Knepley {
576494e7359SMatthew G. Knepley   PetscReal apb, pn1, pn2;
577494e7359SMatthew G. Knepley   PetscInt  k;
578494e7359SMatthew G. Knepley 
579494e7359SMatthew G. Knepley   PetscFunctionBegin;
580494e7359SMatthew G. Knepley   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
581494e7359SMatthew G. Knepley   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
582494e7359SMatthew G. Knepley   apb = a + b;
583494e7359SMatthew G. Knepley   pn2 = 1.0;
584494e7359SMatthew G. Knepley   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
585494e7359SMatthew G. Knepley   *P  = 0.0;
586494e7359SMatthew G. Knepley   for (k = 2; k < n+1; ++k) {
587494e7359SMatthew G. Knepley     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
588494e7359SMatthew G. Knepley     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
589494e7359SMatthew G. Knepley     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
590494e7359SMatthew G. Knepley     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
591494e7359SMatthew G. Knepley 
592494e7359SMatthew G. Knepley     a2  = a2 / a1;
593494e7359SMatthew G. Knepley     a3  = a3 / a1;
594494e7359SMatthew G. Knepley     a4  = a4 / a1;
595494e7359SMatthew G. Knepley     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
596494e7359SMatthew G. Knepley     pn2 = pn1;
597494e7359SMatthew G. Knepley     pn1 = *P;
598494e7359SMatthew G. Knepley   }
599494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
600494e7359SMatthew G. Knepley }
601494e7359SMatthew G. Knepley 
602494e7359SMatthew G. Knepley #undef __FUNCT__
603494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative"
604494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
605494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
606494e7359SMatthew G. Knepley {
607494e7359SMatthew G. Knepley   PetscReal      nP;
608494e7359SMatthew G. Knepley   PetscErrorCode ierr;
609494e7359SMatthew G. Knepley 
610494e7359SMatthew G. Knepley   PetscFunctionBegin;
611494e7359SMatthew G. Knepley   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
612494e7359SMatthew G. Knepley   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
613494e7359SMatthew G. Knepley   *P   = 0.5 * (a + b + n + 1) * nP;
614494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
615494e7359SMatthew G. Knepley }
616494e7359SMatthew G. Knepley 
617494e7359SMatthew G. Knepley #undef __FUNCT__
618494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
619494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
620494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
621494e7359SMatthew G. Knepley {
622494e7359SMatthew G. Knepley   PetscFunctionBegin;
623494e7359SMatthew G. Knepley   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
624494e7359SMatthew G. Knepley   *eta = y;
625494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
626494e7359SMatthew G. Knepley }
627494e7359SMatthew G. Knepley 
628494e7359SMatthew G. Knepley #undef __FUNCT__
629494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
630494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
631494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
632494e7359SMatthew G. Knepley {
633494e7359SMatthew G. Knepley   PetscFunctionBegin;
634494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
635494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
636494e7359SMatthew G. Knepley   *zeta = z;
637494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
638494e7359SMatthew G. Knepley }
639494e7359SMatthew G. Knepley 
640494e7359SMatthew G. Knepley #undef __FUNCT__
641494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
642494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
643494e7359SMatthew G. Knepley {
644494e7359SMatthew G. Knepley   PetscInt       maxIter = 100;
645494e7359SMatthew G. Knepley   PetscReal      eps     = 1.0e-8;
646a8291ba1SSatish Balay   PetscReal      a1, a2, a3, a4, a5, a6;
647494e7359SMatthew G. Knepley   PetscInt       k;
648494e7359SMatthew G. Knepley   PetscErrorCode ierr;
649494e7359SMatthew G. Knepley 
650494e7359SMatthew G. Knepley   PetscFunctionBegin;
651a8291ba1SSatish Balay 
6528b49ba18SBarry Smith   a1      = PetscPowReal(2.0, a+b+1);
653a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA)
6540646a658SBarry Smith   a2      = PetscTGamma(a + npoints + 1);
6550646a658SBarry Smith   a3      = PetscTGamma(b + npoints + 1);
6560646a658SBarry Smith   a4      = PetscTGamma(a + b + npoints + 1);
657a8291ba1SSatish Balay #else
658a8291ba1SSatish Balay   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
659a8291ba1SSatish Balay #endif
660a8291ba1SSatish Balay 
661494e7359SMatthew G. Knepley   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
662494e7359SMatthew G. Knepley   a6   = a1 * a2 * a3 / a4 / a5;
663494e7359SMatthew G. Knepley   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
664494e7359SMatthew G. Knepley    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
665494e7359SMatthew G. Knepley   for (k = 0; k < npoints; ++k) {
6668b49ba18SBarry Smith     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
667494e7359SMatthew G. Knepley     PetscInt  j;
668494e7359SMatthew G. Knepley 
669494e7359SMatthew G. Knepley     if (k > 0) r = 0.5 * (r + x[k-1]);
670494e7359SMatthew G. Knepley     for (j = 0; j < maxIter; ++j) {
671494e7359SMatthew G. Knepley       PetscReal s = 0.0, delta, f, fp;
672494e7359SMatthew G. Knepley       PetscInt  i;
673494e7359SMatthew G. Knepley 
674494e7359SMatthew G. Knepley       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
675494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
676494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
677494e7359SMatthew G. Knepley       delta = f / (fp - f * s);
678494e7359SMatthew G. Knepley       r     = r - delta;
67977b4d14cSPeter Brune       if (PetscAbsReal(delta) < eps) break;
680494e7359SMatthew G. Knepley     }
681494e7359SMatthew G. Knepley     x[k] = r;
682494e7359SMatthew G. Knepley     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
683494e7359SMatthew G. Knepley     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
684494e7359SMatthew G. Knepley   }
685494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
686494e7359SMatthew G. Knepley }
687494e7359SMatthew G. Knepley 
688494e7359SMatthew G. Knepley #undef __FUNCT__
689494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
690fd9d31fbSMatthew G. Knepley /*@C
691494e7359SMatthew G. Knepley   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
692494e7359SMatthew G. Knepley 
693494e7359SMatthew G. Knepley   Not Collective
694494e7359SMatthew G. Knepley 
695494e7359SMatthew G. Knepley   Input Arguments:
696494e7359SMatthew G. Knepley + dim   - The simplex dimension
697744bafbcSMatthew G. Knepley . order - The number of points in one dimension
698494e7359SMatthew G. Knepley . a     - left end of interval (often-1)
699494e7359SMatthew G. Knepley - b     - right end of interval (often +1)
700494e7359SMatthew G. Knepley 
701744bafbcSMatthew G. Knepley   Output Argument:
702552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
703494e7359SMatthew G. Knepley 
704494e7359SMatthew G. Knepley   Level: intermediate
705494e7359SMatthew G. Knepley 
706494e7359SMatthew G. Knepley   References:
707494e7359SMatthew G. Knepley   Karniadakis and Sherwin.
708494e7359SMatthew G. Knepley   FIAT
709494e7359SMatthew G. Knepley 
710744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
711494e7359SMatthew G. Knepley @*/
712552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
713494e7359SMatthew G. Knepley {
714552aa4f7SMatthew G. Knepley   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
715494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
716494e7359SMatthew G. Knepley   PetscInt       i, j, k;
717494e7359SMatthew G. Knepley   PetscErrorCode ierr;
718494e7359SMatthew G. Knepley 
719494e7359SMatthew G. Knepley   PetscFunctionBegin;
720494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
721785e854fSJed Brown   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
722785e854fSJed Brown   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
723494e7359SMatthew G. Knepley   switch (dim) {
724707aa5c5SMatthew G. Knepley   case 0:
725707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
726707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
727785e854fSJed Brown     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
728785e854fSJed Brown     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
729707aa5c5SMatthew G. Knepley     x[0] = 0.0;
730707aa5c5SMatthew G. Knepley     w[0] = 1.0;
731707aa5c5SMatthew G. Knepley     break;
732494e7359SMatthew G. Knepley   case 1:
733552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
734494e7359SMatthew G. Knepley     break;
735494e7359SMatthew G. Knepley   case 2:
736dcca6d9dSJed Brown     ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr);
737552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
738552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
739552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
740552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
741552aa4f7SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
742552aa4f7SMatthew G. Knepley         w[i*order+j] = 0.5 * wx[i] * wy[j];
743494e7359SMatthew G. Knepley       }
744494e7359SMatthew G. Knepley     }
745494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
746494e7359SMatthew G. Knepley     break;
747494e7359SMatthew G. Knepley   case 3:
748dcca6d9dSJed Brown     ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr);
749552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
750552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
751552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
752552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
753552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
754552aa4f7SMatthew G. Knepley         for (k = 0; k < order; ++k) {
755552aa4f7SMatthew 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);
756552aa4f7SMatthew G. Knepley           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
757494e7359SMatthew G. Knepley         }
758494e7359SMatthew G. Knepley       }
759494e7359SMatthew G. Knepley     }
760494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
761494e7359SMatthew G. Knepley     break;
762494e7359SMatthew G. Knepley   default:
763494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
764494e7359SMatthew G. Knepley   }
76521454ff5SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
766bcede257SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*q, order);CHKERRQ(ierr);
76721454ff5SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
768494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
769494e7359SMatthew G. Knepley }
770494e7359SMatthew G. Knepley 
771494e7359SMatthew G. Knepley #undef __FUNCT__
772b3c0f97bSTom Klotz #define __FUNCT__ "PetscDTTanhSinhTensorQuadrature"
773b3c0f97bSTom Klotz /*@C
774b3c0f97bSTom Klotz   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
775b3c0f97bSTom Klotz 
776b3c0f97bSTom Klotz   Not Collective
777b3c0f97bSTom Klotz 
778b3c0f97bSTom Klotz   Input Arguments:
779b3c0f97bSTom Klotz + dim   - The cell dimension
780b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l
781b3c0f97bSTom Klotz . a     - left end of interval (often-1)
782b3c0f97bSTom Klotz - b     - right end of interval (often +1)
783b3c0f97bSTom Klotz 
784b3c0f97bSTom Klotz   Output Argument:
785b3c0f97bSTom Klotz . q - A PetscQuadrature object
786b3c0f97bSTom Klotz 
787b3c0f97bSTom Klotz   Level: intermediate
788b3c0f97bSTom Klotz 
789b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature()
790b3c0f97bSTom Klotz @*/
791b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
792b3c0f97bSTom Klotz {
793b3c0f97bSTom Klotz   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
794b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;                  /* Half-width of the integration interval */
795b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;                  /* Center of the integration interval */
796b3c0f97bSTom Klotz   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
797b3c0f97bSTom Klotz   PetscReal       wk    = 0.5*PETSC_PI;              /* Quadrature weight at x_k */
798b3c0f97bSTom Klotz   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
799b3c0f97bSTom Klotz   PetscReal      *x, *w;
800b3c0f97bSTom Klotz   PetscInt        K, k, npoints;
801b3c0f97bSTom Klotz   PetscErrorCode  ierr;
802b3c0f97bSTom Klotz 
803b3c0f97bSTom Klotz   PetscFunctionBegin;
804b3c0f97bSTom Klotz   if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
805b3c0f97bSTom Klotz   if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
806b3c0f97bSTom Klotz   /* Find K such that the weights are < 32 digits of precision */
807b3c0f97bSTom Klotz   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
808*9add2064SThomas Klotz     wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
809b3c0f97bSTom Klotz   }
810b3c0f97bSTom Klotz   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
811b3c0f97bSTom Klotz   ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr);
812b3c0f97bSTom Klotz   npoints = 2*K-1;
813b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
814b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
815b3c0f97bSTom Klotz   /* Center term */
816b3c0f97bSTom Klotz   x[0] = beta;
817b3c0f97bSTom Klotz   w[0] = 0.5*alpha*PETSC_PI;
818b3c0f97bSTom Klotz   for (k = 1; k < K; ++k) {
819*9add2064SThomas Klotz     wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
820*9add2064SThomas Klotz     xk = tanh(0.5*PETSC_PI*PetscSinhReal(k*h));
821b3c0f97bSTom Klotz     x[2*k-1] = -alpha*xk+beta;
822b3c0f97bSTom Klotz     w[2*k-1] = wk;
823b3c0f97bSTom Klotz     x[2*k+0] =  alpha*xk+beta;
824b3c0f97bSTom Klotz     w[2*k+0] = wk;
825b3c0f97bSTom Klotz   }
826b3c0f97bSTom Klotz   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
827b3c0f97bSTom Klotz   PetscFunctionReturn(0);
828b3c0f97bSTom Klotz }
829b3c0f97bSTom Klotz 
830b3c0f97bSTom Klotz #undef __FUNCT__
831b3c0f97bSTom Klotz #define __FUNCT__ "PetscDTTanhSinhIntegrate"
832b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
833b3c0f97bSTom Klotz {
834b3c0f97bSTom Klotz   const PetscInt  p     = 16;        /* Digits of precision in the evaluation */
835b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;  /* Half-width of the integration interval */
836b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;  /* Center of the integration interval */
837b3c0f97bSTom Klotz   PetscReal       h     = 1.0;       /* Step size, length between x_k */
838b3c0f97bSTom Klotz   PetscInt        l     = 0;         /* Level of refinement, h = 2^{-l} */
839b3c0f97bSTom Klotz   PetscReal       osum  = 0.0;       /* Integral on last level */
840b3c0f97bSTom Klotz   PetscReal       psum  = 0.0;       /* Integral on the level before the last level */
841b3c0f97bSTom Klotz   PetscReal       sum;               /* Integral on current level */
842b3c0f97bSTom Klotz   PetscReal       xk;                /* Quadrature point x_k on reference domain [-1, 1] */
843b3c0f97bSTom Klotz   PetscReal       lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
844b3c0f97bSTom Klotz   PetscReal       wk;                /* Quadrature weight at x_k */
845b3c0f97bSTom Klotz   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
846b3c0f97bSTom Klotz   PetscInt        d;                 /* Digits of precision in the integral */
847b3c0f97bSTom Klotz   PetscErrorCode  ierr;
848b3c0f97bSTom Klotz 
849b3c0f97bSTom Klotz   PetscFunctionBegin;
850b3c0f97bSTom Klotz   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
851b3c0f97bSTom Klotz   /* Center term */
852b3c0f97bSTom Klotz   func(beta, &lval);
853b3c0f97bSTom Klotz   sum = 0.5*alpha*PETSC_PI*lval;
854b3c0f97bSTom Klotz   /* */
855b3c0f97bSTom Klotz   do {
856b3c0f97bSTom Klotz     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
857b3c0f97bSTom Klotz     PetscInt  k = 1;
858b3c0f97bSTom Klotz 
859b3c0f97bSTom Klotz     ++l;
860b3c0f97bSTom Klotz     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
861b3c0f97bSTom Klotz     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
862b3c0f97bSTom Klotz     psum = osum;
863b3c0f97bSTom Klotz     osum = sum;
864b3c0f97bSTom Klotz     h   *= 0.5;
865b3c0f97bSTom Klotz     sum *= 0.5;
866b3c0f97bSTom Klotz     do {
867*9add2064SThomas Klotz       wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
868*9add2064SThomas Klotz       xk = tanh(0.5*PETSC_PI*PetscSinhReal(k*h));
869b3c0f97bSTom Klotz       lx = -alpha*xk+beta;
870b3c0f97bSTom Klotz       rx =  alpha*xk+beta;
871b3c0f97bSTom Klotz       func(lx, &lval);
872b3c0f97bSTom Klotz       func(rx, &rval);
873b3c0f97bSTom Klotz       lterm   = alpha*wk*lval;
874b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
875b3c0f97bSTom Klotz       sum    += lterm;
876b3c0f97bSTom Klotz       rterm   = alpha*wk*rval;
877b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
878b3c0f97bSTom Klotz       sum    += rterm;
879b3c0f97bSTom Klotz       /* if (l == 1) printf("k is %d and sum is %15.15f and wk is %15.15f\n", k, sum, wk); */
880b3c0f97bSTom Klotz       ++k;
881b3c0f97bSTom Klotz       /* Only need to evaluate every other point on refined levels */
882b3c0f97bSTom Klotz       if (l != 1) ++k;
883*9add2064SThomas Klotz     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
884b3c0f97bSTom Klotz 
885b3c0f97bSTom Klotz     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
886b3c0f97bSTom Klotz     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
887b3c0f97bSTom Klotz     d3 = PetscLog10Real(maxTerm) - p;
888b3c0f97bSTom Klotz     d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
889b3c0f97bSTom Klotz     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
890*9add2064SThomas Klotz   } while (d < digits && l < 12);
891b3c0f97bSTom Klotz   *sol = sum;
892b3c0f97bSTom Klotz   PetscFunctionReturn(0);
893b3c0f97bSTom Klotz }
894b3c0f97bSTom Klotz 
895b3c0f97bSTom Klotz #undef __FUNCT__
896194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR"
897194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
898194825f6SJed Brown  * A in column-major format
899194825f6SJed Brown  * Ainv in row-major format
900194825f6SJed Brown  * tau has length m
901194825f6SJed Brown  * worksize must be >= max(1,n)
902194825f6SJed Brown  */
903194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
904194825f6SJed Brown {
905194825f6SJed Brown   PetscErrorCode ierr;
906194825f6SJed Brown   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
907194825f6SJed Brown   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
908194825f6SJed Brown 
909194825f6SJed Brown   PetscFunctionBegin;
910194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
911194825f6SJed Brown   {
912194825f6SJed Brown     PetscInt i,j;
913dcca6d9dSJed Brown     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
914194825f6SJed Brown     for (j=0; j<n; j++) {
915194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
916194825f6SJed Brown     }
917194825f6SJed Brown     mstride = m;
918194825f6SJed Brown   }
919194825f6SJed Brown #else
920194825f6SJed Brown   A = A_in;
921194825f6SJed Brown   Ainv = Ainv_out;
922194825f6SJed Brown #endif
923194825f6SJed Brown 
924194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
925194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
926194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
927194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
928194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
929001a771dSBarry Smith   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
930194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
931194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
932194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
933194825f6SJed Brown 
934194825f6SJed Brown   /* Extract an explicit representation of Q */
935194825f6SJed Brown   Q = Ainv;
936194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
937194825f6SJed Brown   K = N;                        /* full rank */
938001a771dSBarry Smith   PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
939194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
940194825f6SJed Brown 
941194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
942194825f6SJed Brown   Alpha = 1.0;
943194825f6SJed Brown   ldb = lda;
944001a771dSBarry Smith   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
945194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
946194825f6SJed Brown 
947194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
948194825f6SJed Brown   {
949194825f6SJed Brown     PetscInt i;
950194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
951194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
952194825f6SJed Brown   }
953194825f6SJed Brown #endif
954194825f6SJed Brown   PetscFunctionReturn(0);
955194825f6SJed Brown }
956194825f6SJed Brown 
957194825f6SJed Brown #undef __FUNCT__
958194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate"
959194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
960194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
961194825f6SJed Brown {
962194825f6SJed Brown   PetscErrorCode ierr;
963194825f6SJed Brown   PetscReal      *Bv;
964194825f6SJed Brown   PetscInt       i,j;
965194825f6SJed Brown 
966194825f6SJed Brown   PetscFunctionBegin;
967785e854fSJed Brown   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
968194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
969194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
970194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
971194825f6SJed Brown   for (i=0; i<ninterval; i++) {
972194825f6SJed Brown     for (j=0; j<ndegree; j++) {
973194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
974194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
975194825f6SJed Brown     }
976194825f6SJed Brown   }
977194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
978194825f6SJed Brown   PetscFunctionReturn(0);
979194825f6SJed Brown }
980194825f6SJed Brown 
981194825f6SJed Brown #undef __FUNCT__
982194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly"
983194825f6SJed Brown /*@
984194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
985194825f6SJed Brown 
986194825f6SJed Brown    Not Collective
987194825f6SJed Brown 
988194825f6SJed Brown    Input Arguments:
989194825f6SJed Brown +  degree - degree of reconstruction polynomial
990194825f6SJed Brown .  nsource - number of source intervals
991194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
992194825f6SJed Brown .  ntarget - number of target intervals
993194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
994194825f6SJed Brown 
995194825f6SJed Brown    Output Arguments:
996194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
997194825f6SJed Brown 
998194825f6SJed Brown    Level: advanced
999194825f6SJed Brown 
1000194825f6SJed Brown .seealso: PetscDTLegendreEval()
1001194825f6SJed Brown @*/
1002194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1003194825f6SJed Brown {
1004194825f6SJed Brown   PetscErrorCode ierr;
1005194825f6SJed Brown   PetscInt       i,j,k,*bdegrees,worksize;
1006194825f6SJed Brown   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1007194825f6SJed Brown   PetscScalar    *tau,*work;
1008194825f6SJed Brown 
1009194825f6SJed Brown   PetscFunctionBegin;
1010194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
1011194825f6SJed Brown   PetscValidRealPointer(targetx,5);
1012194825f6SJed Brown   PetscValidRealPointer(R,6);
1013194825f6SJed 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);
1014194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
1015194825f6SJed Brown   for (i=0; i<nsource; i++) {
101657622a8eSBarry 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]);
1017194825f6SJed Brown   }
1018194825f6SJed Brown   for (i=0; i<ntarget; i++) {
101957622a8eSBarry 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]);
1020194825f6SJed Brown   }
1021194825f6SJed Brown #endif
1022194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
1023194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1024194825f6SJed Brown   center = (xmin + xmax)/2;
1025194825f6SJed Brown   hscale = (xmax - xmin)/2;
1026194825f6SJed Brown   worksize = nsource;
1027dcca6d9dSJed Brown   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
1028dcca6d9dSJed Brown   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
1029194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1030194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1031194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
1032194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
1033194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1034194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
1035194825f6SJed Brown   for (i=0; i<ntarget; i++) {
1036194825f6SJed Brown     PetscReal rowsum = 0;
1037194825f6SJed Brown     for (j=0; j<nsource; j++) {
1038194825f6SJed Brown       PetscReal sum = 0;
1039194825f6SJed Brown       for (k=0; k<degree+1; k++) {
1040194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1041194825f6SJed Brown       }
1042194825f6SJed Brown       R[i*nsource+j] = sum;
1043194825f6SJed Brown       rowsum += sum;
1044194825f6SJed Brown     }
1045194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1046194825f6SJed Brown   }
1047194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
1048194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
1049194825f6SJed Brown   PetscFunctionReturn(0);
1050194825f6SJed Brown }
1051