xref: /petsc/src/dm/dt/interface/dt.c (revision c74b4a09d4cc2dcefa4b9a7ab96237243e1a985a)
137045ce4SJed Brown /* Discretization tools */
237045ce4SJed Brown 
30c35b76eSJed Brown #include <petscdt.h>            /*I "petscdt.h" I*/
437045ce4SJed Brown #include <petscblaslapack.h>
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
6af0996ceSBarry Smith #include <petsc/private/dtimpl.h>
7665c2dedSJed Brown #include <petscviewer.h>
859804f93SMatthew G. Knepley #include <petscdmplex.h>
959804f93SMatthew G. Knepley #include <petscdmshell.h>
1037045ce4SJed Brown 
1198c04793SMatthew G. Knepley #if defined(PETSC_HAVE_MPFR)
1298c04793SMatthew G. Knepley #include <mpfr.h>
1398c04793SMatthew G. Knepley #endif
1498c04793SMatthew G. Knepley 
150bfcf5a5SMatthew G. Knepley static PetscBool GaussCite       = PETSC_FALSE;
160bfcf5a5SMatthew G. Knepley const char       GaussCitation[] = "@article{GolubWelsch1969,\n"
170bfcf5a5SMatthew G. Knepley                                    "  author  = {Golub and Welsch},\n"
180bfcf5a5SMatthew G. Knepley                                    "  title   = {Calculation of Quadrature Rules},\n"
190bfcf5a5SMatthew G. Knepley                                    "  journal = {Math. Comp.},\n"
200bfcf5a5SMatthew G. Knepley                                    "  volume  = {23},\n"
210bfcf5a5SMatthew G. Knepley                                    "  number  = {106},\n"
220bfcf5a5SMatthew G. Knepley                                    "  pages   = {221--230},\n"
230bfcf5a5SMatthew G. Knepley                                    "  year    = {1969}\n}\n";
240bfcf5a5SMatthew G. Knepley 
2540d8ff71SMatthew G. Knepley /*@
2640d8ff71SMatthew G. Knepley   PetscQuadratureCreate - Create a PetscQuadrature object
2740d8ff71SMatthew G. Knepley 
2840d8ff71SMatthew G. Knepley   Collective on MPI_Comm
2940d8ff71SMatthew G. Knepley 
3040d8ff71SMatthew G. Knepley   Input Parameter:
3140d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object
3240d8ff71SMatthew G. Knepley 
3340d8ff71SMatthew G. Knepley   Output Parameter:
3440d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
3540d8ff71SMatthew G. Knepley 
3640d8ff71SMatthew G. Knepley   Level: beginner
3740d8ff71SMatthew G. Knepley 
3840d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, create
3940d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
4040d8ff71SMatthew G. Knepley @*/
4121454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
4221454ff5SMatthew G. Knepley {
4321454ff5SMatthew G. Knepley   PetscErrorCode ierr;
4421454ff5SMatthew G. Knepley 
4521454ff5SMatthew G. Knepley   PetscFunctionBegin;
4621454ff5SMatthew G. Knepley   PetscValidPointer(q, 2);
47623436dcSMatthew G. Knepley   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
4873107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*q,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
4921454ff5SMatthew G. Knepley   (*q)->dim       = -1;
50a6b92713SMatthew G. Knepley   (*q)->Nc        =  1;
51bcede257SMatthew G. Knepley   (*q)->order     = -1;
5221454ff5SMatthew G. Knepley   (*q)->numPoints = 0;
5321454ff5SMatthew G. Knepley   (*q)->points    = NULL;
5421454ff5SMatthew G. Knepley   (*q)->weights   = NULL;
5521454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
5621454ff5SMatthew G. Knepley }
5721454ff5SMatthew G. Knepley 
58c9638911SMatthew G. Knepley /*@
59c9638911SMatthew G. Knepley   PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object
60c9638911SMatthew G. Knepley 
61c9638911SMatthew G. Knepley   Collective on PetscQuadrature
62c9638911SMatthew G. Knepley 
63c9638911SMatthew G. Knepley   Input Parameter:
64c9638911SMatthew G. Knepley . q  - The PetscQuadrature object
65c9638911SMatthew G. Knepley 
66c9638911SMatthew G. Knepley   Output Parameter:
67c9638911SMatthew G. Knepley . r  - The new PetscQuadrature object
68c9638911SMatthew G. Knepley 
69c9638911SMatthew G. Knepley   Level: beginner
70c9638911SMatthew G. Knepley 
71c9638911SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, clone
72c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData()
73c9638911SMatthew G. Knepley @*/
74c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
75c9638911SMatthew G. Knepley {
76a6b92713SMatthew G. Knepley   PetscInt         order, dim, Nc, Nq;
77c9638911SMatthew G. Knepley   const PetscReal *points, *weights;
78c9638911SMatthew G. Knepley   PetscReal       *p, *w;
79c9638911SMatthew G. Knepley   PetscErrorCode   ierr;
80c9638911SMatthew G. Knepley 
81c9638911SMatthew G. Knepley   PetscFunctionBegin;
82c9638911SMatthew G. Knepley   PetscValidPointer(q, 2);
83c9638911SMatthew G. Knepley   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr);
84c9638911SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
85c9638911SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr);
86a6b92713SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
87c9638911SMatthew G. Knepley   ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr);
88f0a0bfafSMatthew G. Knepley   ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr);
89c9638911SMatthew G. Knepley   ierr = PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));CHKERRQ(ierr);
90a6b92713SMatthew G. Knepley   ierr = PetscMemcpy(w, weights, Nc * Nq * sizeof(PetscReal));CHKERRQ(ierr);
91a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr);
92c9638911SMatthew G. Knepley   PetscFunctionReturn(0);
93c9638911SMatthew G. Knepley }
94c9638911SMatthew G. Knepley 
9540d8ff71SMatthew G. Knepley /*@
9640d8ff71SMatthew G. Knepley   PetscQuadratureDestroy - Destroys a PetscQuadrature object
9740d8ff71SMatthew G. Knepley 
9840d8ff71SMatthew G. Knepley   Collective on PetscQuadrature
9940d8ff71SMatthew G. Knepley 
10040d8ff71SMatthew G. Knepley   Input Parameter:
10140d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
10240d8ff71SMatthew G. Knepley 
10340d8ff71SMatthew G. Knepley   Level: beginner
10440d8ff71SMatthew G. Knepley 
10540d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, destroy
10640d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
10740d8ff71SMatthew G. Knepley @*/
108bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
109bfa639d9SMatthew G. Knepley {
110bfa639d9SMatthew G. Knepley   PetscErrorCode ierr;
111bfa639d9SMatthew G. Knepley 
112bfa639d9SMatthew G. Knepley   PetscFunctionBegin;
11321454ff5SMatthew G. Knepley   if (!*q) PetscFunctionReturn(0);
11421454ff5SMatthew G. Knepley   PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1);
11521454ff5SMatthew G. Knepley   if (--((PetscObject)(*q))->refct > 0) {
11621454ff5SMatthew G. Knepley     *q = NULL;
11721454ff5SMatthew G. Knepley     PetscFunctionReturn(0);
11821454ff5SMatthew G. Knepley   }
11921454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
12021454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
12121454ff5SMatthew G. Knepley   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
12221454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
12321454ff5SMatthew G. Knepley }
12421454ff5SMatthew G. Knepley 
125bcede257SMatthew G. Knepley /*@
126a6b92713SMatthew G. Knepley   PetscQuadratureGetOrder - Return the order of the method
127bcede257SMatthew G. Knepley 
128bcede257SMatthew G. Knepley   Not collective
129bcede257SMatthew G. Knepley 
130bcede257SMatthew G. Knepley   Input Parameter:
131bcede257SMatthew G. Knepley . q - The PetscQuadrature object
132bcede257SMatthew G. Knepley 
133bcede257SMatthew G. Knepley   Output Parameter:
134bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
135bcede257SMatthew G. Knepley 
136bcede257SMatthew G. Knepley   Level: intermediate
137bcede257SMatthew G. Knepley 
138bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
139bcede257SMatthew G. Knepley @*/
140bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
141bcede257SMatthew G. Knepley {
142bcede257SMatthew G. Knepley   PetscFunctionBegin;
143bcede257SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
144bcede257SMatthew G. Knepley   PetscValidPointer(order, 2);
145bcede257SMatthew G. Knepley   *order = q->order;
146bcede257SMatthew G. Knepley   PetscFunctionReturn(0);
147bcede257SMatthew G. Knepley }
148bcede257SMatthew G. Knepley 
149bcede257SMatthew G. Knepley /*@
150a6b92713SMatthew G. Knepley   PetscQuadratureSetOrder - Return the order of the method
151bcede257SMatthew G. Knepley 
152bcede257SMatthew G. Knepley   Not collective
153bcede257SMatthew G. Knepley 
154bcede257SMatthew G. Knepley   Input Parameters:
155bcede257SMatthew G. Knepley + q - The PetscQuadrature object
156bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
157bcede257SMatthew G. Knepley 
158bcede257SMatthew G. Knepley   Level: intermediate
159bcede257SMatthew G. Knepley 
160bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
161bcede257SMatthew G. Knepley @*/
162bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
163bcede257SMatthew G. Knepley {
164bcede257SMatthew G. Knepley   PetscFunctionBegin;
165bcede257SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
166bcede257SMatthew G. Knepley   q->order = order;
167bcede257SMatthew G. Knepley   PetscFunctionReturn(0);
168bcede257SMatthew G. Knepley }
169bcede257SMatthew G. Knepley 
170a6b92713SMatthew G. Knepley /*@
171a6b92713SMatthew G. Knepley   PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated
172a6b92713SMatthew G. Knepley 
173a6b92713SMatthew G. Knepley   Not collective
174a6b92713SMatthew G. Knepley 
175a6b92713SMatthew G. Knepley   Input Parameter:
176a6b92713SMatthew G. Knepley . q - The PetscQuadrature object
177a6b92713SMatthew G. Knepley 
178a6b92713SMatthew G. Knepley   Output Parameter:
179a6b92713SMatthew G. Knepley . Nc - The number of components
180a6b92713SMatthew G. Knepley 
181a6b92713SMatthew G. Knepley   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
182a6b92713SMatthew G. Knepley 
183a6b92713SMatthew G. Knepley   Level: intermediate
184a6b92713SMatthew G. Knepley 
185a6b92713SMatthew G. Knepley .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
186a6b92713SMatthew G. Knepley @*/
187a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc)
188a6b92713SMatthew G. Knepley {
189a6b92713SMatthew G. Knepley   PetscFunctionBegin;
190a6b92713SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
191a6b92713SMatthew G. Knepley   PetscValidPointer(Nc, 2);
192a6b92713SMatthew G. Knepley   *Nc = q->Nc;
193a6b92713SMatthew G. Knepley   PetscFunctionReturn(0);
194a6b92713SMatthew G. Knepley }
195a6b92713SMatthew G. Knepley 
196a6b92713SMatthew G. Knepley /*@
197a6b92713SMatthew G. Knepley   PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated
198a6b92713SMatthew G. Knepley 
199a6b92713SMatthew G. Knepley   Not collective
200a6b92713SMatthew G. Knepley 
201a6b92713SMatthew G. Knepley   Input Parameters:
202a6b92713SMatthew G. Knepley + q  - The PetscQuadrature object
203a6b92713SMatthew G. Knepley - Nc - The number of components
204a6b92713SMatthew G. Knepley 
205a6b92713SMatthew G. Knepley   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
206a6b92713SMatthew G. Knepley 
207a6b92713SMatthew G. Knepley   Level: intermediate
208a6b92713SMatthew G. Knepley 
209a6b92713SMatthew G. Knepley .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
210a6b92713SMatthew G. Knepley @*/
211a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc)
212a6b92713SMatthew G. Knepley {
213a6b92713SMatthew G. Knepley   PetscFunctionBegin;
214a6b92713SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
215a6b92713SMatthew G. Knepley   q->Nc = Nc;
216a6b92713SMatthew G. Knepley   PetscFunctionReturn(0);
217a6b92713SMatthew G. Knepley }
218a6b92713SMatthew G. Knepley 
21940d8ff71SMatthew G. Knepley /*@C
22040d8ff71SMatthew G. Knepley   PetscQuadratureGetData - Returns the data defining the quadrature
22140d8ff71SMatthew G. Knepley 
22240d8ff71SMatthew G. Knepley   Not collective
22340d8ff71SMatthew G. Knepley 
22440d8ff71SMatthew G. Knepley   Input Parameter:
22540d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
22640d8ff71SMatthew G. Knepley 
22740d8ff71SMatthew G. Knepley   Output Parameters:
22840d8ff71SMatthew G. Knepley + dim - The spatial dimension
229805e7170SToby Isaac . Nc - The number of components
23040d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
23140d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
23240d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
23340d8ff71SMatthew G. Knepley 
23440d8ff71SMatthew G. Knepley   Level: intermediate
23540d8ff71SMatthew G. Knepley 
23695452b02SPatrick Sanan   Fortran Notes:
23795452b02SPatrick Sanan     From Fortran you must call PetscQuadratureRestoreData() when you are done with the data
2381fd49c25SBarry Smith 
23940d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature
24040d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
24140d8ff71SMatthew G. Knepley @*/
242a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
24321454ff5SMatthew G. Knepley {
24421454ff5SMatthew G. Knepley   PetscFunctionBegin;
24521454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
24621454ff5SMatthew G. Knepley   if (dim) {
24721454ff5SMatthew G. Knepley     PetscValidPointer(dim, 2);
24821454ff5SMatthew G. Knepley     *dim = q->dim;
24921454ff5SMatthew G. Knepley   }
250a6b92713SMatthew G. Knepley   if (Nc) {
251a6b92713SMatthew G. Knepley     PetscValidPointer(Nc, 3);
252a6b92713SMatthew G. Knepley     *Nc = q->Nc;
253a6b92713SMatthew G. Knepley   }
25421454ff5SMatthew G. Knepley   if (npoints) {
255a6b92713SMatthew G. Knepley     PetscValidPointer(npoints, 4);
25621454ff5SMatthew G. Knepley     *npoints = q->numPoints;
25721454ff5SMatthew G. Knepley   }
25821454ff5SMatthew G. Knepley   if (points) {
259a6b92713SMatthew G. Knepley     PetscValidPointer(points, 5);
26021454ff5SMatthew G. Knepley     *points = q->points;
26121454ff5SMatthew G. Knepley   }
26221454ff5SMatthew G. Knepley   if (weights) {
263a6b92713SMatthew G. Knepley     PetscValidPointer(weights, 6);
26421454ff5SMatthew G. Knepley     *weights = q->weights;
26521454ff5SMatthew G. Knepley   }
26621454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
26721454ff5SMatthew G. Knepley }
26821454ff5SMatthew G. Knepley 
26940d8ff71SMatthew G. Knepley /*@C
27040d8ff71SMatthew G. Knepley   PetscQuadratureSetData - Sets the data defining the quadrature
27140d8ff71SMatthew G. Knepley 
27240d8ff71SMatthew G. Knepley   Not collective
27340d8ff71SMatthew G. Knepley 
27440d8ff71SMatthew G. Knepley   Input Parameters:
27540d8ff71SMatthew G. Knepley + q  - The PetscQuadrature object
27640d8ff71SMatthew G. Knepley . dim - The spatial dimension
277e2b35d93SBarry Smith . Nc - The number of components
27840d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
27940d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
28040d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
28140d8ff71SMatthew G. Knepley 
282c99e0549SMatthew G. Knepley   Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them.
283f2fd9e53SMatthew G. Knepley 
28440d8ff71SMatthew G. Knepley   Level: intermediate
28540d8ff71SMatthew G. Knepley 
28640d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature
28740d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
28840d8ff71SMatthew G. Knepley @*/
289a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
29021454ff5SMatthew G. Knepley {
29121454ff5SMatthew G. Knepley   PetscFunctionBegin;
29221454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
29321454ff5SMatthew G. Knepley   if (dim >= 0)     q->dim       = dim;
294a6b92713SMatthew G. Knepley   if (Nc >= 0)      q->Nc        = Nc;
29521454ff5SMatthew G. Knepley   if (npoints >= 0) q->numPoints = npoints;
29621454ff5SMatthew G. Knepley   if (points) {
29721454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
29821454ff5SMatthew G. Knepley     q->points = points;
29921454ff5SMatthew G. Knepley   }
30021454ff5SMatthew G. Knepley   if (weights) {
30121454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
30221454ff5SMatthew G. Knepley     q->weights = weights;
30321454ff5SMatthew G. Knepley   }
304f9fd7fdbSMatthew G. Knepley   PetscFunctionReturn(0);
305f9fd7fdbSMatthew G. Knepley }
306f9fd7fdbSMatthew G. Knepley 
307d9bac1caSLisandro Dalcin static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v)
308d9bac1caSLisandro Dalcin {
309d9bac1caSLisandro Dalcin   PetscInt          q, d, c;
310d9bac1caSLisandro Dalcin   PetscViewerFormat format;
311d9bac1caSLisandro Dalcin   PetscErrorCode    ierr;
312d9bac1caSLisandro Dalcin 
313d9bac1caSLisandro Dalcin   PetscFunctionBegin;
314*c74b4a09SMatthew G. Knepley   if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D) with %D components\n", quad->order, quad->numPoints, quad->dim, quad->Nc);CHKERRQ(ierr);}
315*c74b4a09SMatthew G. Knepley   else              {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);}
316d9bac1caSLisandro Dalcin   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
317d9bac1caSLisandro Dalcin   if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
318d9bac1caSLisandro Dalcin   for (q = 0; q < quad->numPoints; ++q) {
319*c74b4a09SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr);
320d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr);
321d9bac1caSLisandro Dalcin     for (d = 0; d < quad->dim; ++d) {
322d9bac1caSLisandro Dalcin       if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
323d9bac1caSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
324d9bac1caSLisandro Dalcin     }
325d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr);
326*c74b4a09SMatthew G. Knepley     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);}
327d9bac1caSLisandro Dalcin     for (c = 0; c < quad->Nc; ++c) {
328d9bac1caSLisandro Dalcin       if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
329*c74b4a09SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr);
330d9bac1caSLisandro Dalcin     }
331d9bac1caSLisandro Dalcin     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);}
332d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr);
333d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr);
334d9bac1caSLisandro Dalcin   }
335d9bac1caSLisandro Dalcin   PetscFunctionReturn(0);
336d9bac1caSLisandro Dalcin }
337d9bac1caSLisandro Dalcin 
33840d8ff71SMatthew G. Knepley /*@C
33940d8ff71SMatthew G. Knepley   PetscQuadratureView - Views a PetscQuadrature object
34040d8ff71SMatthew G. Knepley 
34140d8ff71SMatthew G. Knepley   Collective on PetscQuadrature
34240d8ff71SMatthew G. Knepley 
34340d8ff71SMatthew G. Knepley   Input Parameters:
344d9bac1caSLisandro Dalcin + quad  - The PetscQuadrature object
34540d8ff71SMatthew G. Knepley - viewer - The PetscViewer object
34640d8ff71SMatthew G. Knepley 
34740d8ff71SMatthew G. Knepley   Level: beginner
34840d8ff71SMatthew G. Knepley 
34940d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, view
35040d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
35140d8ff71SMatthew G. Knepley @*/
352f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
353f9fd7fdbSMatthew G. Knepley {
354d9bac1caSLisandro Dalcin   PetscBool      iascii;
355f9fd7fdbSMatthew G. Knepley   PetscErrorCode ierr;
356f9fd7fdbSMatthew G. Knepley 
357f9fd7fdbSMatthew G. Knepley   PetscFunctionBegin;
358d9bac1caSLisandro Dalcin   PetscValidHeader(quad, 1);
359d9bac1caSLisandro Dalcin   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
360d9bac1caSLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);}
361d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
362d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
363d9bac1caSLisandro Dalcin   if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);}
364d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
365bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
366bfa639d9SMatthew G. Knepley }
367bfa639d9SMatthew G. Knepley 
36889710940SMatthew G. Knepley /*@C
36989710940SMatthew G. Knepley   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
37089710940SMatthew G. Knepley 
37189710940SMatthew G. Knepley   Not collective
37289710940SMatthew G. Knepley 
37389710940SMatthew G. Knepley   Input Parameter:
37489710940SMatthew G. Knepley + q - The original PetscQuadrature
37589710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into
37689710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement
37789710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement
37889710940SMatthew G. Knepley 
37989710940SMatthew G. Knepley   Output Parameters:
38089710940SMatthew G. Knepley . dim - The dimension
38189710940SMatthew G. Knepley 
38289710940SMatthew G. Knepley   Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
38389710940SMatthew G. Knepley 
384f5f57ec0SBarry Smith  Not available from Fortran
385f5f57ec0SBarry Smith 
38689710940SMatthew G. Knepley   Level: intermediate
38789710940SMatthew G. Knepley 
38889710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
38989710940SMatthew G. Knepley @*/
39089710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
39189710940SMatthew G. Knepley {
39289710940SMatthew G. Knepley   const PetscReal *points,    *weights;
39389710940SMatthew G. Knepley   PetscReal       *pointsRef, *weightsRef;
394a6b92713SMatthew G. Knepley   PetscInt         dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
39589710940SMatthew G. Knepley   PetscErrorCode   ierr;
39689710940SMatthew G. Knepley 
39789710940SMatthew G. Knepley   PetscFunctionBegin;
39889710940SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
39989710940SMatthew G. Knepley   PetscValidPointer(v0, 3);
40089710940SMatthew G. Knepley   PetscValidPointer(jac, 4);
40189710940SMatthew G. Knepley   PetscValidPointer(qref, 5);
40289710940SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr);
40389710940SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
404a6b92713SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr);
40589710940SMatthew G. Knepley   npointsRef = npoints*numSubelements;
40689710940SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr);
407a6b92713SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr);
40889710940SMatthew G. Knepley   for (c = 0; c < numSubelements; ++c) {
40989710940SMatthew G. Knepley     for (p = 0; p < npoints; ++p) {
41089710940SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
41189710940SMatthew G. Knepley         pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
41289710940SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
41389710940SMatthew G. Knepley           pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
41489710940SMatthew G. Knepley         }
41589710940SMatthew G. Knepley       }
41689710940SMatthew G. Knepley       /* Could also use detJ here */
417a6b92713SMatthew G. Knepley       for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements;
41889710940SMatthew G. Knepley     }
41989710940SMatthew G. Knepley   }
42089710940SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr);
421a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr);
42289710940SMatthew G. Knepley   PetscFunctionReturn(0);
42389710940SMatthew G. Knepley }
42489710940SMatthew G. Knepley 
42537045ce4SJed Brown /*@
42637045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
42737045ce4SJed Brown 
42837045ce4SJed Brown    Not Collective
42937045ce4SJed Brown 
43037045ce4SJed Brown    Input Arguments:
43137045ce4SJed Brown +  npoints - number of spatial points to evaluate at
43237045ce4SJed Brown .  points - array of locations to evaluate at
43337045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
43437045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
43537045ce4SJed Brown 
43637045ce4SJed Brown    Output Arguments:
4370298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
4380298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
4390298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
44037045ce4SJed Brown 
44137045ce4SJed Brown    Level: intermediate
44237045ce4SJed Brown 
44337045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
44437045ce4SJed Brown @*/
44537045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
44637045ce4SJed Brown {
44737045ce4SJed Brown   PetscInt i,maxdegree;
44837045ce4SJed Brown 
44937045ce4SJed Brown   PetscFunctionBegin;
45037045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
45137045ce4SJed Brown   maxdegree = degrees[ndegree-1];
45237045ce4SJed Brown   for (i=0; i<npoints; i++) {
45337045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
45437045ce4SJed Brown     PetscInt  j,k;
45537045ce4SJed Brown     x    = points[i];
45637045ce4SJed Brown     pm2  = 0;
45737045ce4SJed Brown     pm1  = 1;
45837045ce4SJed Brown     pd2  = 0;
45937045ce4SJed Brown     pd1  = 0;
46037045ce4SJed Brown     pdd2 = 0;
46137045ce4SJed Brown     pdd1 = 0;
46237045ce4SJed Brown     k    = 0;
46337045ce4SJed Brown     if (degrees[k] == 0) {
46437045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
46537045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
46637045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
46737045ce4SJed Brown       k++;
46837045ce4SJed Brown     }
46937045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
47037045ce4SJed Brown       PetscReal p,d,dd;
47137045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
47237045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
47337045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
47437045ce4SJed Brown       pm2  = pm1;
47537045ce4SJed Brown       pm1  = p;
47637045ce4SJed Brown       pd2  = pd1;
47737045ce4SJed Brown       pd1  = d;
47837045ce4SJed Brown       pdd2 = pdd1;
47937045ce4SJed Brown       pdd1 = dd;
48037045ce4SJed Brown       if (degrees[k] == j) {
48137045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
48237045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
48337045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
48437045ce4SJed Brown       }
48537045ce4SJed Brown     }
48637045ce4SJed Brown   }
48737045ce4SJed Brown   PetscFunctionReturn(0);
48837045ce4SJed Brown }
48937045ce4SJed Brown 
49037045ce4SJed Brown /*@
49137045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
49237045ce4SJed Brown 
49337045ce4SJed Brown    Not Collective
49437045ce4SJed Brown 
49537045ce4SJed Brown    Input Arguments:
49637045ce4SJed Brown +  npoints - number of points
49737045ce4SJed Brown .  a - left end of interval (often-1)
49837045ce4SJed Brown -  b - right end of interval (often +1)
49937045ce4SJed Brown 
50037045ce4SJed Brown    Output Arguments:
50137045ce4SJed Brown +  x - quadrature points
50237045ce4SJed Brown -  w - quadrature weights
50337045ce4SJed Brown 
50437045ce4SJed Brown    Level: intermediate
50537045ce4SJed Brown 
50637045ce4SJed Brown    References:
50796a0c994SBarry Smith .   1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
50837045ce4SJed Brown 
50937045ce4SJed Brown .seealso: PetscDTLegendreEval()
51037045ce4SJed Brown @*/
51137045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
51237045ce4SJed Brown {
51337045ce4SJed Brown   PetscErrorCode ierr;
51437045ce4SJed Brown   PetscInt       i;
51537045ce4SJed Brown   PetscReal      *work;
51637045ce4SJed Brown   PetscScalar    *Z;
51737045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
51837045ce4SJed Brown 
51937045ce4SJed Brown   PetscFunctionBegin;
5200bfcf5a5SMatthew G. Knepley   ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr);
52137045ce4SJed Brown   /* Set up the Golub-Welsch system */
52237045ce4SJed Brown   for (i=0; i<npoints; i++) {
52337045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
52437045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
52537045ce4SJed Brown   }
526dcca6d9dSJed Brown   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
527c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
52837045ce4SJed Brown   LDZ  = N;
52937045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
5308b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
53137045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
5321c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
53337045ce4SJed Brown 
53437045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
53537045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
53637045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
53719a57d60SBarry Smith     if (x[i] == -0.0) x[i] = 0.0;
53837045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
5390d644c17SKarl Rupp 
54088393a60SJed Brown     w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
54137045ce4SJed Brown   }
54237045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
54337045ce4SJed Brown   PetscFunctionReturn(0);
54437045ce4SJed Brown }
545194825f6SJed Brown 
5468272889dSSatish Balay static void qAndLEvaluation(PetscInt n, PetscReal x, PetscReal *q, PetscReal *qp, PetscReal *Ln)
5478272889dSSatish Balay /*
5488272889dSSatish Balay   Compute the polynomial q(x) = L_{N+1}(x) - L_{n-1}(x) and its derivative in
5498272889dSSatish Balay   addition to L_N(x) as these are needed for computing the GLL points via Newton's method.
5508272889dSSatish Balay   Reference: "Implementing Spectral Methods for Partial Differential Equations: Algorithms
5518272889dSSatish Balay   for Scientists and Engineers" by David A. Kopriva.
5528272889dSSatish Balay */
5538272889dSSatish Balay {
5548272889dSSatish Balay   PetscInt k;
5558272889dSSatish Balay 
5568272889dSSatish Balay   PetscReal Lnp;
5578272889dSSatish Balay   PetscReal Lnp1, Lnp1p;
5588272889dSSatish Balay   PetscReal Lnm1, Lnm1p;
5598272889dSSatish Balay   PetscReal Lnm2, Lnm2p;
5608272889dSSatish Balay 
5618272889dSSatish Balay   Lnm1  = 1.0;
5628272889dSSatish Balay   *Ln   = x;
5638272889dSSatish Balay   Lnm1p = 0.0;
5648272889dSSatish Balay   Lnp   = 1.0;
5658272889dSSatish Balay 
5668272889dSSatish Balay   for (k=2; k<=n; ++k) {
5678272889dSSatish Balay     Lnm2  = Lnm1;
5688272889dSSatish Balay     Lnm1  = *Ln;
5698272889dSSatish Balay     Lnm2p = Lnm1p;
5708272889dSSatish Balay     Lnm1p = Lnp;
5718272889dSSatish Balay     *Ln   = (2.*((PetscReal)k)-1.)/(1.0*((PetscReal)k))*x*Lnm1 - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm2;
5728272889dSSatish Balay     Lnp   = Lnm2p + (2.0*((PetscReal)k)-1.)*Lnm1;
5738272889dSSatish Balay   }
5748272889dSSatish Balay   k     = n+1;
5758272889dSSatish Balay   Lnp1  = (2.*((PetscReal)k)-1.)/(((PetscReal)k))*x*(*Ln) - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm1;
5768272889dSSatish Balay   Lnp1p = Lnm1p + (2.0*((PetscReal)k)-1.)*(*Ln);
5778272889dSSatish Balay   *q    = Lnp1 - Lnm1;
5788272889dSSatish Balay   *qp   = Lnp1p - Lnm1p;
5798272889dSSatish Balay }
5808272889dSSatish Balay 
5818272889dSSatish Balay /*@C
5828272889dSSatish Balay    PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre
5838272889dSSatish Balay                       nodes of a given size on the domain [-1,1]
5848272889dSSatish Balay 
5858272889dSSatish Balay    Not Collective
5868272889dSSatish Balay 
5878272889dSSatish Balay    Input Parameter:
5888272889dSSatish Balay +  n - number of grid nodes
589f2e8fe4dShannah_mairs -  type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
5908272889dSSatish Balay 
5918272889dSSatish Balay    Output Arguments:
5928272889dSSatish Balay +  x - quadrature points
5938272889dSSatish Balay -  w - quadrature weights
5948272889dSSatish Balay 
5958272889dSSatish Balay    Notes:
5968272889dSSatish Balay     For n > 30  the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
5978272889dSSatish Balay           close enough to the desired solution
5988272889dSSatish Balay 
5998272889dSSatish Balay    These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes
6008272889dSSatish Balay 
6018272889dSSatish Balay    See  http://epubs.siam.org/doi/abs/10.1137/110855442  http://epubs.siam.org/doi/abs/10.1137/120889873 for better ways to compute GLL nodes
6028272889dSSatish Balay 
6038272889dSSatish Balay    Level: intermediate
6048272889dSSatish Balay 
6058272889dSSatish Balay .seealso: PetscDTGaussQuadrature()
6068272889dSSatish Balay 
6078272889dSSatish Balay @*/
608916e780bShannah_mairs PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w)
6098272889dSSatish Balay {
6108272889dSSatish Balay   PetscErrorCode ierr;
6118272889dSSatish Balay 
6128272889dSSatish Balay   PetscFunctionBegin;
6138272889dSSatish Balay   if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element");
6148272889dSSatish Balay 
615f2e8fe4dShannah_mairs   if (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA) {
6168272889dSSatish Balay     PetscReal      *M,si;
6178272889dSSatish Balay     PetscBLASInt   bn,lierr;
6188272889dSSatish Balay     PetscReal      x0,z0,z1,z2;
6198272889dSSatish Balay     PetscInt       i,p = npoints - 1,nn;
6208272889dSSatish Balay 
6218272889dSSatish Balay     x[0]   =-1.0;
6228272889dSSatish Balay     x[npoints-1] = 1.0;
6238272889dSSatish Balay     if (npoints-2 > 0){
6248272889dSSatish Balay       ierr = PetscMalloc1(npoints-1,&M);CHKERRQ(ierr);
6258272889dSSatish Balay       for (i=0; i<npoints-2; i++) {
6268272889dSSatish Balay         si  = ((PetscReal)i)+1.0;
6278272889dSSatish Balay         M[i]=0.5*PetscSqrtReal(si*(si+2.0)/((si+0.5)*(si+1.5)));
6288272889dSSatish Balay       }
6298272889dSSatish Balay       ierr = PetscBLASIntCast(npoints-2,&bn);CHKERRQ(ierr);
6308272889dSSatish Balay       ierr = PetscMemzero(&x[1],bn*sizeof(x[1]));CHKERRQ(ierr);
6318272889dSSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
6328272889dSSatish Balay       x0=0;
6338272889dSSatish Balay       PetscStackCallBLAS("LAPACKsteqr",LAPACKREALsteqr_("N",&bn,&x[1],M,&x0,&bn,M,&lierr));
6348272889dSSatish Balay       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in STERF Lapack routine %d",(int)lierr);
6358272889dSSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
6368272889dSSatish Balay       ierr = PetscFree(M);CHKERRQ(ierr);
6378272889dSSatish Balay     }
6388272889dSSatish Balay     if ((npoints-1)%2==0) {
6398272889dSSatish Balay       x[(npoints-1)/2]   = 0.0; /* hard wire to exactly 0.0 since linear algebra produces nonzero */
6408272889dSSatish Balay     }
6418272889dSSatish Balay 
6428272889dSSatish Balay     w[0] = w[p] = 2.0/(((PetscReal)(p))*(((PetscReal)p)+1.0));
6438272889dSSatish Balay     z2 = -1.;                      /* Dummy value to avoid -Wmaybe-initialized */
6448272889dSSatish Balay     for (i=1; i<p; i++) {
6458272889dSSatish Balay       x0  = x[i];
6468272889dSSatish Balay       z0 = 1.0;
6478272889dSSatish Balay       z1 = x0;
6488272889dSSatish Balay       for (nn=1; nn<p; nn++) {
6498272889dSSatish Balay         z2 = x0*z1*(2.0*((PetscReal)nn)+1.0)/(((PetscReal)nn)+1.0)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.0));
6508272889dSSatish Balay         z0 = z1;
6518272889dSSatish Balay         z1 = z2;
6528272889dSSatish Balay       }
6538272889dSSatish Balay       w[i]=2.0/(((PetscReal)p)*(((PetscReal)p)+1.0)*z2*z2);
6548272889dSSatish Balay     }
6558272889dSSatish Balay   } else {
6568272889dSSatish Balay     PetscInt  j,m;
6578272889dSSatish Balay     PetscReal z1,z,q,qp,Ln;
6588272889dSSatish Balay     PetscReal *pt;
6598272889dSSatish Balay     ierr = PetscMalloc1(npoints,&pt);CHKERRQ(ierr);
6608272889dSSatish Balay 
661d410ae54Shannah_mairs     if (npoints > 30) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON produces incorrect answers for n > 30");
6628272889dSSatish Balay     x[0]     = -1.0;
6638272889dSSatish Balay     x[npoints-1]   = 1.0;
6648272889dSSatish Balay     w[0]   = w[npoints-1] = 2./(((PetscReal)npoints)*(((PetscReal)npoints)-1.0));;
6658272889dSSatish Balay     m  = (npoints-1)/2; /* The roots are symmetric, so we only find half of them. */
6668272889dSSatish Balay     for (j=1; j<=m; j++) { /* Loop over the desired roots. */
6678272889dSSatish Balay       z = -1.0*PetscCosReal((PETSC_PI*((PetscReal)j)+0.25)/(((PetscReal)npoints)-1.0))-(3.0/(8.0*(((PetscReal)npoints)-1.0)*PETSC_PI))*(1.0/(((PetscReal)j)+0.25));
6688272889dSSatish Balay       /* Starting with the above approximation to the ith root, we enter */
6698272889dSSatish Balay       /* the main loop of refinement by Newton's method.                 */
6708272889dSSatish Balay       do {
6718272889dSSatish Balay         qAndLEvaluation(npoints-1,z,&q,&qp,&Ln);
6728272889dSSatish Balay         z1 = z;
6738272889dSSatish Balay         z  = z1-q/qp; /* Newton's method. */
6748272889dSSatish Balay       } while (PetscAbs(z-z1) > 10.*PETSC_MACHINE_EPSILON);
6758272889dSSatish Balay       qAndLEvaluation(npoints-1,z,&q,&qp,&Ln);
6768272889dSSatish Balay 
6778272889dSSatish Balay       x[j]       = z;
6788272889dSSatish Balay       x[npoints-1-j]   = -z;      /* and put in its symmetric counterpart.   */
6798272889dSSatish Balay       w[j]     = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln);  /* Compute the weight */
6808272889dSSatish Balay       w[npoints-1-j] = w[j];                 /* and its symmetric counterpart. */
6818272889dSSatish Balay       pt[j]=qp;
6828272889dSSatish Balay     }
6838272889dSSatish Balay 
6848272889dSSatish Balay     if ((npoints-1)%2==0) {
6858272889dSSatish Balay       qAndLEvaluation(npoints-1,0.0,&q,&qp,&Ln);
6868272889dSSatish Balay       x[(npoints-1)/2]   = 0.0;
6878272889dSSatish Balay       w[(npoints-1)/2] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln);
6888272889dSSatish Balay     }
6898272889dSSatish Balay     ierr = PetscFree(pt);CHKERRQ(ierr);
6908272889dSSatish Balay   }
6918272889dSSatish Balay   PetscFunctionReturn(0);
6928272889dSSatish Balay }
6938272889dSSatish Balay 
694744bafbcSMatthew G. Knepley /*@
695744bafbcSMatthew G. Knepley   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
696744bafbcSMatthew G. Knepley 
697744bafbcSMatthew G. Knepley   Not Collective
698744bafbcSMatthew G. Knepley 
699744bafbcSMatthew G. Knepley   Input Arguments:
700744bafbcSMatthew G. Knepley + dim     - The spatial dimension
701a6b92713SMatthew G. Knepley . Nc      - The number of components
702744bafbcSMatthew G. Knepley . npoints - number of points in one dimension
703744bafbcSMatthew G. Knepley . a       - left end of interval (often-1)
704744bafbcSMatthew G. Knepley - b       - right end of interval (often +1)
705744bafbcSMatthew G. Knepley 
706744bafbcSMatthew G. Knepley   Output Argument:
707744bafbcSMatthew G. Knepley . q - A PetscQuadrature object
708744bafbcSMatthew G. Knepley 
709744bafbcSMatthew G. Knepley   Level: intermediate
710744bafbcSMatthew G. Knepley 
711744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
712744bafbcSMatthew G. Knepley @*/
713a6b92713SMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
714744bafbcSMatthew G. Knepley {
715a6b92713SMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
716744bafbcSMatthew G. Knepley   PetscReal     *x, *w, *xw, *ww;
717744bafbcSMatthew G. Knepley   PetscErrorCode ierr;
718744bafbcSMatthew G. Knepley 
719744bafbcSMatthew G. Knepley   PetscFunctionBegin;
720744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
721a6b92713SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr);
722744bafbcSMatthew G. Knepley   /* Set up the Golub-Welsch system */
723744bafbcSMatthew G. Knepley   switch (dim) {
724744bafbcSMatthew G. Knepley   case 0:
725744bafbcSMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
726744bafbcSMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
727744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
728a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
729744bafbcSMatthew G. Knepley     x[0] = 0.0;
730a6b92713SMatthew G. Knepley     for (c = 0; c < Nc; ++c) w[c] = 1.0;
731744bafbcSMatthew G. Knepley     break;
732744bafbcSMatthew G. Knepley   case 1:
733a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr);
734a6b92713SMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr);
735a6b92713SMatthew G. Knepley     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i];
736a6b92713SMatthew G. Knepley     ierr = PetscFree(ww);CHKERRQ(ierr);
737744bafbcSMatthew G. Knepley     break;
738744bafbcSMatthew G. Knepley   case 2:
739744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
740744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
741744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
742744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
743744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+0] = xw[i];
744744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+1] = xw[j];
745a6b92713SMatthew G. Knepley         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j];
746744bafbcSMatthew G. Knepley       }
747744bafbcSMatthew G. Knepley     }
748744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
749744bafbcSMatthew G. Knepley     break;
750744bafbcSMatthew G. Knepley   case 3:
751744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
752744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
753744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
754744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
755744bafbcSMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
756744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
757744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
758744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
759a6b92713SMatthew G. Knepley           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k];
760744bafbcSMatthew G. Knepley         }
761744bafbcSMatthew G. Knepley       }
762744bafbcSMatthew G. Knepley     }
763744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
764744bafbcSMatthew G. Knepley     break;
765744bafbcSMatthew G. Knepley   default:
766744bafbcSMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
767744bafbcSMatthew G. Knepley   }
768744bafbcSMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
7692f5fb066SToby Isaac   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
770a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
771d9bac1caSLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr);
772744bafbcSMatthew G. Knepley   PetscFunctionReturn(0);
773744bafbcSMatthew G. Knepley }
774744bafbcSMatthew G. Knepley 
775494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
776494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
777494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
778494e7359SMatthew G. Knepley {
779494e7359SMatthew G. Knepley   PetscReal f = 1.0;
780494e7359SMatthew G. Knepley   PetscInt  i;
781494e7359SMatthew G. Knepley 
782494e7359SMatthew G. Knepley   PetscFunctionBegin;
783494e7359SMatthew G. Knepley   for (i = 1; i < n+1; ++i) f *= i;
784494e7359SMatthew G. Knepley   *factorial = f;
785494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
786494e7359SMatthew G. Knepley }
787494e7359SMatthew G. Knepley 
788494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
789494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
790494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
791494e7359SMatthew G. Knepley {
792494e7359SMatthew G. Knepley   PetscReal apb, pn1, pn2;
793494e7359SMatthew G. Knepley   PetscInt  k;
794494e7359SMatthew G. Knepley 
795494e7359SMatthew G. Knepley   PetscFunctionBegin;
796494e7359SMatthew G. Knepley   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
797494e7359SMatthew G. Knepley   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
798494e7359SMatthew G. Knepley   apb = a + b;
799494e7359SMatthew G. Knepley   pn2 = 1.0;
800494e7359SMatthew G. Knepley   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
801494e7359SMatthew G. Knepley   *P  = 0.0;
802494e7359SMatthew G. Knepley   for (k = 2; k < n+1; ++k) {
803494e7359SMatthew G. Knepley     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
804494e7359SMatthew G. Knepley     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
805494e7359SMatthew G. Knepley     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
806494e7359SMatthew G. Knepley     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
807494e7359SMatthew G. Knepley 
808494e7359SMatthew G. Knepley     a2  = a2 / a1;
809494e7359SMatthew G. Knepley     a3  = a3 / a1;
810494e7359SMatthew G. Knepley     a4  = a4 / a1;
811494e7359SMatthew G. Knepley     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
812494e7359SMatthew G. Knepley     pn2 = pn1;
813494e7359SMatthew G. Knepley     pn1 = *P;
814494e7359SMatthew G. Knepley   }
815494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
816494e7359SMatthew G. Knepley }
817494e7359SMatthew G. Knepley 
818494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
819494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
820494e7359SMatthew G. Knepley {
821494e7359SMatthew G. Knepley   PetscReal      nP;
822494e7359SMatthew G. Knepley   PetscErrorCode ierr;
823494e7359SMatthew G. Knepley 
824494e7359SMatthew G. Knepley   PetscFunctionBegin;
825494e7359SMatthew G. Knepley   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
826494e7359SMatthew G. Knepley   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
827494e7359SMatthew G. Knepley   *P   = 0.5 * (a + b + n + 1) * nP;
828494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
829494e7359SMatthew G. Knepley }
830494e7359SMatthew G. Knepley 
831494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
832494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
833494e7359SMatthew G. Knepley {
834494e7359SMatthew G. Knepley   PetscFunctionBegin;
835494e7359SMatthew G. Knepley   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
836494e7359SMatthew G. Knepley   *eta = y;
837494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
838494e7359SMatthew G. Knepley }
839494e7359SMatthew G. Knepley 
840494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
841494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
842494e7359SMatthew G. Knepley {
843494e7359SMatthew G. Knepley   PetscFunctionBegin;
844494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
845494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
846494e7359SMatthew G. Knepley   *zeta = z;
847494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
848494e7359SMatthew G. Knepley }
849494e7359SMatthew G. Knepley 
850494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
851494e7359SMatthew G. Knepley {
852494e7359SMatthew G. Knepley   PetscInt       maxIter = 100;
853494e7359SMatthew G. Knepley   PetscReal      eps     = 1.0e-8;
854a8291ba1SSatish Balay   PetscReal      a1, a2, a3, a4, a5, a6;
855494e7359SMatthew G. Knepley   PetscInt       k;
856494e7359SMatthew G. Knepley   PetscErrorCode ierr;
857494e7359SMatthew G. Knepley 
858494e7359SMatthew G. Knepley   PetscFunctionBegin;
859a8291ba1SSatish Balay 
8608b49ba18SBarry Smith   a1      = PetscPowReal(2.0, a+b+1);
861a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA)
8620646a658SBarry Smith   a2      = PetscTGamma(a + npoints + 1);
8630646a658SBarry Smith   a3      = PetscTGamma(b + npoints + 1);
8640646a658SBarry Smith   a4      = PetscTGamma(a + b + npoints + 1);
865a8291ba1SSatish Balay #else
86629bcbfd0SToby Isaac   {
867d24bbb91SToby Isaac     PetscInt ia, ib;
86829bcbfd0SToby Isaac 
869d24bbb91SToby Isaac     ia = (PetscInt) a;
870d24bbb91SToby Isaac     ib = (PetscInt) b;
871d24bbb91SToby 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 */
872d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr);
873d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr);
874d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr);
87529bcbfd0SToby Isaac     } else {
876a8291ba1SSatish Balay       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
87729bcbfd0SToby Isaac     }
87829bcbfd0SToby Isaac   }
879a8291ba1SSatish Balay #endif
880a8291ba1SSatish Balay 
881494e7359SMatthew G. Knepley   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
882494e7359SMatthew G. Knepley   a6   = a1 * a2 * a3 / a4 / a5;
883494e7359SMatthew G. Knepley   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
884494e7359SMatthew G. Knepley    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
885494e7359SMatthew G. Knepley   for (k = 0; k < npoints; ++k) {
8868b49ba18SBarry Smith     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
887494e7359SMatthew G. Knepley     PetscInt  j;
888494e7359SMatthew G. Knepley 
889494e7359SMatthew G. Knepley     if (k > 0) r = 0.5 * (r + x[k-1]);
890494e7359SMatthew G. Knepley     for (j = 0; j < maxIter; ++j) {
891494e7359SMatthew G. Knepley       PetscReal s = 0.0, delta, f, fp;
892494e7359SMatthew G. Knepley       PetscInt  i;
893494e7359SMatthew G. Knepley 
894494e7359SMatthew G. Knepley       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
895494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
896494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
897494e7359SMatthew G. Knepley       delta = f / (fp - f * s);
898494e7359SMatthew G. Knepley       r     = r - delta;
89977b4d14cSPeter Brune       if (PetscAbsReal(delta) < eps) break;
900494e7359SMatthew G. Knepley     }
901494e7359SMatthew G. Knepley     x[k] = r;
902494e7359SMatthew G. Knepley     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
903494e7359SMatthew G. Knepley     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
904494e7359SMatthew G. Knepley   }
905494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
906494e7359SMatthew G. Knepley }
907494e7359SMatthew G. Knepley 
908f5f57ec0SBarry Smith /*@
909494e7359SMatthew G. Knepley   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
910494e7359SMatthew G. Knepley 
911494e7359SMatthew G. Knepley   Not Collective
912494e7359SMatthew G. Knepley 
913494e7359SMatthew G. Knepley   Input Arguments:
914494e7359SMatthew G. Knepley + dim     - The simplex dimension
915a6b92713SMatthew G. Knepley . Nc      - The number of components
916dcce0ee2SMatthew G. Knepley . npoints - The number of points in one dimension
917494e7359SMatthew G. Knepley . a       - left end of interval (often-1)
918494e7359SMatthew G. Knepley - b       - right end of interval (often +1)
919494e7359SMatthew G. Knepley 
920744bafbcSMatthew G. Knepley   Output Argument:
921552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
922494e7359SMatthew G. Knepley 
923494e7359SMatthew G. Knepley   Level: intermediate
924494e7359SMatthew G. Knepley 
925494e7359SMatthew G. Knepley   References:
92696a0c994SBarry Smith .  1. - Karniadakis and Sherwin.  FIAT
927494e7359SMatthew G. Knepley 
928744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
929494e7359SMatthew G. Knepley @*/
930dcce0ee2SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
931494e7359SMatthew G. Knepley {
932dcce0ee2SMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints;
933494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
934a6b92713SMatthew G. Knepley   PetscInt       i, j, k, c;
935494e7359SMatthew G. Knepley   PetscErrorCode ierr;
936494e7359SMatthew G. Knepley 
937494e7359SMatthew G. Knepley   PetscFunctionBegin;
938494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
939dcce0ee2SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr);
940dcce0ee2SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr);
941494e7359SMatthew G. Knepley   switch (dim) {
942707aa5c5SMatthew G. Knepley   case 0:
943707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
944707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
945785e854fSJed Brown     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
946a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
947707aa5c5SMatthew G. Knepley     x[0] = 0.0;
948a6b92713SMatthew G. Knepley     for (c = 0; c < Nc; ++c) w[c] = 1.0;
949707aa5c5SMatthew G. Knepley     break;
950494e7359SMatthew G. Knepley   case 1:
951dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr);
952dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr);
953dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i];
954a6b92713SMatthew G. Knepley     ierr = PetscFree(wx);CHKERRQ(ierr);
955494e7359SMatthew G. Knepley     break;
956494e7359SMatthew G. Knepley   case 2:
957dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr);
958dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
959dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
960dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
961dcce0ee2SMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
962dcce0ee2SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr);
963dcce0ee2SMatthew G. Knepley         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j];
964494e7359SMatthew G. Knepley       }
965494e7359SMatthew G. Knepley     }
966494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
967494e7359SMatthew G. Knepley     break;
968494e7359SMatthew G. Knepley   case 3:
969dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr);
970dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
971dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
972dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
973dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
974dcce0ee2SMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
975dcce0ee2SMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
976dcce0ee2SMatthew G. Knepley           ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*npoints+j)*npoints+k)*3+0], &x[((i*npoints+j)*npoints+k)*3+1], &x[((i*npoints+j)*npoints+k)*3+2]);CHKERRQ(ierr);
977dcce0ee2SMatthew G. Knepley           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k];
978494e7359SMatthew G. Knepley         }
979494e7359SMatthew G. Knepley       }
980494e7359SMatthew G. Knepley     }
981494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
982494e7359SMatthew G. Knepley     break;
983494e7359SMatthew G. Knepley   default:
984494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
985494e7359SMatthew G. Knepley   }
98621454ff5SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
9872f5fb066SToby Isaac   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
988dcce0ee2SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
989d9bac1caSLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr);
990494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
991494e7359SMatthew G. Knepley }
992494e7359SMatthew G. Knepley 
993f5f57ec0SBarry Smith /*@
994b3c0f97bSTom Klotz   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
995b3c0f97bSTom Klotz 
996b3c0f97bSTom Klotz   Not Collective
997b3c0f97bSTom Klotz 
998b3c0f97bSTom Klotz   Input Arguments:
999b3c0f97bSTom Klotz + dim   - The cell dimension
1000b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l
1001b3c0f97bSTom Klotz . a     - left end of interval (often-1)
1002b3c0f97bSTom Klotz - b     - right end of interval (often +1)
1003b3c0f97bSTom Klotz 
1004b3c0f97bSTom Klotz   Output Argument:
1005b3c0f97bSTom Klotz . q - A PetscQuadrature object
1006b3c0f97bSTom Klotz 
1007b3c0f97bSTom Klotz   Level: intermediate
1008b3c0f97bSTom Klotz 
1009b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature()
1010b3c0f97bSTom Klotz @*/
1011b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
1012b3c0f97bSTom Klotz {
1013b3c0f97bSTom Klotz   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
1014b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;                  /* Half-width of the integration interval */
1015b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;                  /* Center of the integration interval */
1016b3c0f97bSTom Klotz   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
1017d84b4d08SMatthew G. Knepley   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
1018b3c0f97bSTom Klotz   PetscReal       wk    = 0.5*PETSC_PI;              /* Quadrature weight at x_k */
1019b3c0f97bSTom Klotz   PetscReal      *x, *w;
1020b3c0f97bSTom Klotz   PetscInt        K, k, npoints;
1021b3c0f97bSTom Klotz   PetscErrorCode  ierr;
1022b3c0f97bSTom Klotz 
1023b3c0f97bSTom Klotz   PetscFunctionBegin;
1024b3c0f97bSTom Klotz   if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
1025b3c0f97bSTom Klotz   if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
1026b3c0f97bSTom Klotz   /* Find K such that the weights are < 32 digits of precision */
1027b3c0f97bSTom Klotz   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
10289add2064SThomas Klotz     wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
1029b3c0f97bSTom Klotz   }
1030b3c0f97bSTom Klotz   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
1031b3c0f97bSTom Klotz   ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr);
1032b3c0f97bSTom Klotz   npoints = 2*K-1;
1033b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
1034b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
1035b3c0f97bSTom Klotz   /* Center term */
1036b3c0f97bSTom Klotz   x[0] = beta;
1037b3c0f97bSTom Klotz   w[0] = 0.5*alpha*PETSC_PI;
1038b3c0f97bSTom Klotz   for (k = 1; k < K; ++k) {
10399add2064SThomas Klotz     wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
10401118d4bcSLisandro Dalcin     xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h));
1041b3c0f97bSTom Klotz     x[2*k-1] = -alpha*xk+beta;
1042b3c0f97bSTom Klotz     w[2*k-1] = wk;
1043b3c0f97bSTom Klotz     x[2*k+0] =  alpha*xk+beta;
1044b3c0f97bSTom Klotz     w[2*k+0] = wk;
1045b3c0f97bSTom Klotz   }
1046a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr);
1047b3c0f97bSTom Klotz   PetscFunctionReturn(0);
1048b3c0f97bSTom Klotz }
1049b3c0f97bSTom Klotz 
1050b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1051b3c0f97bSTom Klotz {
1052b3c0f97bSTom Klotz   const PetscInt  p     = 16;        /* Digits of precision in the evaluation */
1053b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;  /* Half-width of the integration interval */
1054b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;  /* Center of the integration interval */
1055b3c0f97bSTom Klotz   PetscReal       h     = 1.0;       /* Step size, length between x_k */
1056b3c0f97bSTom Klotz   PetscInt        l     = 0;         /* Level of refinement, h = 2^{-l} */
1057b3c0f97bSTom Klotz   PetscReal       osum  = 0.0;       /* Integral on last level */
1058b3c0f97bSTom Klotz   PetscReal       psum  = 0.0;       /* Integral on the level before the last level */
1059b3c0f97bSTom Klotz   PetscReal       sum;               /* Integral on current level */
1060446c295cSMatthew G. Knepley   PetscReal       yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1061b3c0f97bSTom Klotz   PetscReal       lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1062b3c0f97bSTom Klotz   PetscReal       wk;                /* Quadrature weight at x_k */
1063b3c0f97bSTom Klotz   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
1064b3c0f97bSTom Klotz   PetscInt        d;                 /* Digits of precision in the integral */
1065b3c0f97bSTom Klotz 
1066b3c0f97bSTom Klotz   PetscFunctionBegin;
1067b3c0f97bSTom Klotz   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1068b3c0f97bSTom Klotz   /* Center term */
1069b3c0f97bSTom Klotz   func(beta, &lval);
1070b3c0f97bSTom Klotz   sum = 0.5*alpha*PETSC_PI*lval;
1071b3c0f97bSTom Klotz   /* */
1072b3c0f97bSTom Klotz   do {
1073b3c0f97bSTom Klotz     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
1074b3c0f97bSTom Klotz     PetscInt  k = 1;
1075b3c0f97bSTom Klotz 
1076b3c0f97bSTom Klotz     ++l;
1077b3c0f97bSTom Klotz     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1078b3c0f97bSTom Klotz     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1079b3c0f97bSTom Klotz     psum = osum;
1080b3c0f97bSTom Klotz     osum = sum;
1081b3c0f97bSTom Klotz     h   *= 0.5;
1082b3c0f97bSTom Klotz     sum *= 0.5;
1083b3c0f97bSTom Klotz     do {
10849add2064SThomas Klotz       wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1085446c295cSMatthew G. Knepley       yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1086446c295cSMatthew G. Knepley       lx = -alpha*(1.0 - yk)+beta;
1087446c295cSMatthew G. Knepley       rx =  alpha*(1.0 - yk)+beta;
1088b3c0f97bSTom Klotz       func(lx, &lval);
1089b3c0f97bSTom Klotz       func(rx, &rval);
1090b3c0f97bSTom Klotz       lterm   = alpha*wk*lval;
1091b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
1092b3c0f97bSTom Klotz       sum    += lterm;
1093b3c0f97bSTom Klotz       rterm   = alpha*wk*rval;
1094b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
1095b3c0f97bSTom Klotz       sum    += rterm;
1096b3c0f97bSTom Klotz       ++k;
1097b3c0f97bSTom Klotz       /* Only need to evaluate every other point on refined levels */
1098b3c0f97bSTom Klotz       if (l != 1) ++k;
10999add2064SThomas Klotz     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
1100b3c0f97bSTom Klotz 
1101b3c0f97bSTom Klotz     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
1102b3c0f97bSTom Klotz     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
1103b3c0f97bSTom Klotz     d3 = PetscLog10Real(maxTerm) - p;
110409d48545SBarry Smith     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
110509d48545SBarry Smith     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
1106b3c0f97bSTom Klotz     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
11079add2064SThomas Klotz   } while (d < digits && l < 12);
1108b3c0f97bSTom Klotz   *sol = sum;
1109e510cb1fSThomas Klotz 
1110b3c0f97bSTom Klotz   PetscFunctionReturn(0);
1111b3c0f97bSTom Klotz }
1112b3c0f97bSTom Klotz 
1113497880caSRichard Tran Mills #if defined(PETSC_HAVE_MPFR)
111429f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
111529f144ccSMatthew G. Knepley {
1116e510cb1fSThomas Klotz   const PetscInt  safetyFactor = 2;  /* Calculate abcissa until 2*p digits */
111729f144ccSMatthew G. Knepley   PetscInt        l            = 0;  /* Level of refinement, h = 2^{-l} */
111829f144ccSMatthew G. Knepley   mpfr_t          alpha;             /* Half-width of the integration interval */
111929f144ccSMatthew G. Knepley   mpfr_t          beta;              /* Center of the integration interval */
112029f144ccSMatthew G. Knepley   mpfr_t          h;                 /* Step size, length between x_k */
112129f144ccSMatthew G. Knepley   mpfr_t          osum;              /* Integral on last level */
112229f144ccSMatthew G. Knepley   mpfr_t          psum;              /* Integral on the level before the last level */
112329f144ccSMatthew G. Knepley   mpfr_t          sum;               /* Integral on current level */
112429f144ccSMatthew G. Knepley   mpfr_t          yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
112529f144ccSMatthew G. Knepley   mpfr_t          lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
112629f144ccSMatthew G. Knepley   mpfr_t          wk;                /* Quadrature weight at x_k */
112729f144ccSMatthew G. Knepley   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
112829f144ccSMatthew G. Knepley   PetscInt        d;                 /* Digits of precision in the integral */
112929f144ccSMatthew G. Knepley   mpfr_t          pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
113029f144ccSMatthew G. Knepley 
113129f144ccSMatthew G. Knepley   PetscFunctionBegin;
113229f144ccSMatthew G. Knepley   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
113329f144ccSMatthew G. Knepley   /* Create high precision storage */
1134c9f744b5SMatthew 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);
113529f144ccSMatthew G. Knepley   /* Initialization */
113629f144ccSMatthew G. Knepley   mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
113729f144ccSMatthew G. Knepley   mpfr_set_d(beta,  0.5*(b+a), MPFR_RNDN);
113829f144ccSMatthew G. Knepley   mpfr_set_d(osum,  0.0,       MPFR_RNDN);
113929f144ccSMatthew G. Knepley   mpfr_set_d(psum,  0.0,       MPFR_RNDN);
114029f144ccSMatthew G. Knepley   mpfr_set_d(h,     1.0,       MPFR_RNDN);
114129f144ccSMatthew G. Knepley   mpfr_const_pi(pi2, MPFR_RNDN);
114229f144ccSMatthew G. Knepley   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
114329f144ccSMatthew G. Knepley   /* Center term */
114429f144ccSMatthew G. Knepley   func(0.5*(b+a), &lval);
114529f144ccSMatthew G. Knepley   mpfr_set(sum, pi2, MPFR_RNDN);
114629f144ccSMatthew G. Knepley   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
114729f144ccSMatthew G. Knepley   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
114829f144ccSMatthew G. Knepley   /* */
114929f144ccSMatthew G. Knepley   do {
115029f144ccSMatthew G. Knepley     PetscReal d1, d2, d3, d4;
115129f144ccSMatthew G. Knepley     PetscInt  k = 1;
115229f144ccSMatthew G. Knepley 
115329f144ccSMatthew G. Knepley     ++l;
115429f144ccSMatthew G. Knepley     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
115529f144ccSMatthew G. Knepley     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
115629f144ccSMatthew G. Knepley     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
115729f144ccSMatthew G. Knepley     mpfr_set(psum, osum, MPFR_RNDN);
115829f144ccSMatthew G. Knepley     mpfr_set(osum,  sum, MPFR_RNDN);
115929f144ccSMatthew G. Knepley     mpfr_mul_d(h,   h,   0.5, MPFR_RNDN);
116029f144ccSMatthew G. Knepley     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
116129f144ccSMatthew G. Knepley     do {
116229f144ccSMatthew G. Knepley       mpfr_set_si(kh, k, MPFR_RNDN);
116329f144ccSMatthew G. Knepley       mpfr_mul(kh, kh, h, MPFR_RNDN);
116429f144ccSMatthew G. Knepley       /* Weight */
116529f144ccSMatthew G. Knepley       mpfr_set(wk, h, MPFR_RNDN);
116629f144ccSMatthew G. Knepley       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
116729f144ccSMatthew G. Knepley       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
116829f144ccSMatthew G. Knepley       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
116929f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
117029f144ccSMatthew G. Knepley       mpfr_sqr(tmp, tmp, MPFR_RNDN);
117129f144ccSMatthew G. Knepley       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
117229f144ccSMatthew G. Knepley       mpfr_div(wk, wk, tmp, MPFR_RNDN);
117329f144ccSMatthew G. Knepley       /* Abscissa */
117429f144ccSMatthew G. Knepley       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
117529f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
117629f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
117729f144ccSMatthew G. Knepley       mpfr_exp(tmp, msinh, MPFR_RNDN);
117829f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
117929f144ccSMatthew G. Knepley       /* Quadrature points */
118029f144ccSMatthew G. Knepley       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
118129f144ccSMatthew G. Knepley       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
118229f144ccSMatthew G. Knepley       mpfr_add(lx, lx, beta, MPFR_RNDU);
118329f144ccSMatthew G. Knepley       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
118429f144ccSMatthew G. Knepley       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
118529f144ccSMatthew G. Knepley       mpfr_add(rx, rx, beta, MPFR_RNDD);
118629f144ccSMatthew G. Knepley       /* Evaluation */
118729f144ccSMatthew G. Knepley       func(mpfr_get_d(lx, MPFR_RNDU), &lval);
118829f144ccSMatthew G. Knepley       func(mpfr_get_d(rx, MPFR_RNDD), &rval);
118929f144ccSMatthew G. Knepley       /* Update */
119029f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
119129f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
119229f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
119329f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
119429f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
119529f144ccSMatthew G. Knepley       mpfr_set(curTerm, tmp, MPFR_RNDN);
119629f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
119729f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
119829f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
119929f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
120029f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
120129f144ccSMatthew G. Knepley       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
120229f144ccSMatthew G. Knepley       ++k;
120329f144ccSMatthew G. Knepley       /* Only need to evaluate every other point on refined levels */
120429f144ccSMatthew G. Knepley       if (l != 1) ++k;
120529f144ccSMatthew G. Knepley       mpfr_log10(tmp, wk, MPFR_RNDN);
120629f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
1207c9f744b5SMatthew G. Knepley     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
120829f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
120929f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
121029f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
121129f144ccSMatthew G. Knepley     d1 = mpfr_get_d(tmp, MPFR_RNDN);
121229f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
121329f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
121429f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
121529f144ccSMatthew G. Knepley     d2 = mpfr_get_d(tmp, MPFR_RNDN);
121629f144ccSMatthew G. Knepley     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
1217c9f744b5SMatthew G. Knepley     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
121829f144ccSMatthew G. Knepley     mpfr_log10(tmp, curTerm, MPFR_RNDN);
121929f144ccSMatthew G. Knepley     d4 = mpfr_get_d(tmp, MPFR_RNDN);
122029f144ccSMatthew G. Knepley     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1221b0649871SThomas Klotz   } while (d < digits && l < 8);
122229f144ccSMatthew G. Knepley   *sol = mpfr_get_d(sum, MPFR_RNDN);
122329f144ccSMatthew G. Knepley   /* Cleanup */
122429f144ccSMatthew G. Knepley   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
122529f144ccSMatthew G. Knepley   PetscFunctionReturn(0);
122629f144ccSMatthew G. Knepley }
1227d525116cSMatthew G. Knepley #else
1228fbfcfee5SBarry Smith 
1229d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1230d525116cSMatthew G. Knepley {
1231d525116cSMatthew G. Knepley   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
1232d525116cSMatthew G. Knepley }
123329f144ccSMatthew G. Knepley #endif
123429f144ccSMatthew G. Knepley 
1235194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
1236194825f6SJed Brown  * A in column-major format
1237194825f6SJed Brown  * Ainv in row-major format
1238194825f6SJed Brown  * tau has length m
1239194825f6SJed Brown  * worksize must be >= max(1,n)
1240194825f6SJed Brown  */
1241194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1242194825f6SJed Brown {
1243194825f6SJed Brown   PetscErrorCode ierr;
1244194825f6SJed Brown   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
1245194825f6SJed Brown   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
1246194825f6SJed Brown 
1247194825f6SJed Brown   PetscFunctionBegin;
1248194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
1249194825f6SJed Brown   {
1250194825f6SJed Brown     PetscInt i,j;
1251dcca6d9dSJed Brown     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
1252194825f6SJed Brown     for (j=0; j<n; j++) {
1253194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1254194825f6SJed Brown     }
1255194825f6SJed Brown     mstride = m;
1256194825f6SJed Brown   }
1257194825f6SJed Brown #else
1258194825f6SJed Brown   A = A_in;
1259194825f6SJed Brown   Ainv = Ainv_out;
1260194825f6SJed Brown #endif
1261194825f6SJed Brown 
1262194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
1263194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
1264194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
1265194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
1266194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1267001a771dSBarry Smith   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1268194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1269194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1270194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1271194825f6SJed Brown 
1272194825f6SJed Brown   /* Extract an explicit representation of Q */
1273194825f6SJed Brown   Q = Ainv;
1274194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
1275194825f6SJed Brown   K = N;                        /* full rank */
1276c964aadfSJose E. Roman   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1277194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1278194825f6SJed Brown 
1279194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1280194825f6SJed Brown   Alpha = 1.0;
1281194825f6SJed Brown   ldb = lda;
1282001a771dSBarry Smith   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1283194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
1284194825f6SJed Brown 
1285194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
1286194825f6SJed Brown   {
1287194825f6SJed Brown     PetscInt i;
1288194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1289194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
1290194825f6SJed Brown   }
1291194825f6SJed Brown #endif
1292194825f6SJed Brown   PetscFunctionReturn(0);
1293194825f6SJed Brown }
1294194825f6SJed Brown 
1295194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1296194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1297194825f6SJed Brown {
1298194825f6SJed Brown   PetscErrorCode ierr;
1299194825f6SJed Brown   PetscReal      *Bv;
1300194825f6SJed Brown   PetscInt       i,j;
1301194825f6SJed Brown 
1302194825f6SJed Brown   PetscFunctionBegin;
1303785e854fSJed Brown   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
1304194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
1305194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
1306194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1307194825f6SJed Brown   for (i=0; i<ninterval; i++) {
1308194825f6SJed Brown     for (j=0; j<ndegree; j++) {
1309194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1310194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1311194825f6SJed Brown     }
1312194825f6SJed Brown   }
1313194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
1314194825f6SJed Brown   PetscFunctionReturn(0);
1315194825f6SJed Brown }
1316194825f6SJed Brown 
1317194825f6SJed Brown /*@
1318194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
1319194825f6SJed Brown 
1320194825f6SJed Brown    Not Collective
1321194825f6SJed Brown 
1322194825f6SJed Brown    Input Arguments:
1323194825f6SJed Brown +  degree - degree of reconstruction polynomial
1324194825f6SJed Brown .  nsource - number of source intervals
1325194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1326194825f6SJed Brown .  ntarget - number of target intervals
1327194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
1328194825f6SJed Brown 
1329194825f6SJed Brown    Output Arguments:
1330194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
1331194825f6SJed Brown 
1332194825f6SJed Brown    Level: advanced
1333194825f6SJed Brown 
1334194825f6SJed Brown .seealso: PetscDTLegendreEval()
1335194825f6SJed Brown @*/
1336194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1337194825f6SJed Brown {
1338194825f6SJed Brown   PetscErrorCode ierr;
1339194825f6SJed Brown   PetscInt       i,j,k,*bdegrees,worksize;
1340194825f6SJed Brown   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1341194825f6SJed Brown   PetscScalar    *tau,*work;
1342194825f6SJed Brown 
1343194825f6SJed Brown   PetscFunctionBegin;
1344194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
1345194825f6SJed Brown   PetscValidRealPointer(targetx,5);
1346194825f6SJed Brown   PetscValidRealPointer(R,6);
1347194825f6SJed 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);
1348194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
1349194825f6SJed Brown   for (i=0; i<nsource; i++) {
135057622a8eSBarry 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]);
1351194825f6SJed Brown   }
1352194825f6SJed Brown   for (i=0; i<ntarget; i++) {
135357622a8eSBarry 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]);
1354194825f6SJed Brown   }
1355194825f6SJed Brown #endif
1356194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
1357194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1358194825f6SJed Brown   center = (xmin + xmax)/2;
1359194825f6SJed Brown   hscale = (xmax - xmin)/2;
1360194825f6SJed Brown   worksize = nsource;
1361dcca6d9dSJed Brown   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
1362dcca6d9dSJed Brown   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
1363194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1364194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1365194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
1366194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
1367194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1368194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
1369194825f6SJed Brown   for (i=0; i<ntarget; i++) {
1370194825f6SJed Brown     PetscReal rowsum = 0;
1371194825f6SJed Brown     for (j=0; j<nsource; j++) {
1372194825f6SJed Brown       PetscReal sum = 0;
1373194825f6SJed Brown       for (k=0; k<degree+1; k++) {
1374194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1375194825f6SJed Brown       }
1376194825f6SJed Brown       R[i*nsource+j] = sum;
1377194825f6SJed Brown       rowsum += sum;
1378194825f6SJed Brown     }
1379194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1380194825f6SJed Brown   }
1381194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
1382194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
1383194825f6SJed Brown   PetscFunctionReturn(0);
1384194825f6SJed Brown }
1385916e780bShannah_mairs 
1386916e780bShannah_mairs /*@C
1387916e780bShannah_mairs    PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points
1388916e780bShannah_mairs 
1389916e780bShannah_mairs    Not Collective
1390916e780bShannah_mairs 
1391916e780bShannah_mairs    Input Parameter:
1392916e780bShannah_mairs +  n - the number of GLL nodes
1393916e780bShannah_mairs .  nodes - the GLL nodes
1394916e780bShannah_mairs .  weights - the GLL weights
1395916e780bShannah_mairs .  f - the function values at the nodes
1396916e780bShannah_mairs 
1397916e780bShannah_mairs    Output Parameter:
1398916e780bShannah_mairs .  in - the value of the integral
1399916e780bShannah_mairs 
1400916e780bShannah_mairs    Level: beginner
1401916e780bShannah_mairs 
1402916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature()
1403916e780bShannah_mairs 
1404916e780bShannah_mairs @*/
1405916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in)
1406916e780bShannah_mairs {
1407916e780bShannah_mairs   PetscInt          i;
1408916e780bShannah_mairs 
1409916e780bShannah_mairs   PetscFunctionBegin;
1410916e780bShannah_mairs   *in = 0.;
1411916e780bShannah_mairs   for (i=0; i<n; i++) {
1412916e780bShannah_mairs     *in += f[i]*f[i]*weights[i];
1413916e780bShannah_mairs   }
1414916e780bShannah_mairs   PetscFunctionReturn(0);
1415916e780bShannah_mairs }
1416916e780bShannah_mairs 
1417916e780bShannah_mairs /*@C
1418916e780bShannah_mairs    PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element
1419916e780bShannah_mairs 
1420916e780bShannah_mairs    Not Collective
1421916e780bShannah_mairs 
1422916e780bShannah_mairs    Input Parameter:
1423916e780bShannah_mairs +  n - the number of GLL nodes
1424916e780bShannah_mairs .  nodes - the GLL nodes
1425916e780bShannah_mairs .  weights - the GLL weights
1426916e780bShannah_mairs 
1427916e780bShannah_mairs    Output Parameter:
1428916e780bShannah_mairs .  A - the stiffness element
1429916e780bShannah_mairs 
1430916e780bShannah_mairs    Level: beginner
1431916e780bShannah_mairs 
1432916e780bShannah_mairs    Notes:
1433916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy()
1434916e780bShannah_mairs 
1435916e780bShannah_mairs    You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented (the array is symmetric)
1436916e780bShannah_mairs 
1437916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1438916e780bShannah_mairs 
1439916e780bShannah_mairs @*/
1440916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1441916e780bShannah_mairs {
1442916e780bShannah_mairs   PetscReal        **A;
1443916e780bShannah_mairs   PetscErrorCode  ierr;
1444916e780bShannah_mairs   const PetscReal  *gllnodes = nodes;
1445916e780bShannah_mairs   const PetscInt   p = n-1;
1446916e780bShannah_mairs   PetscReal        z0,z1,z2 = -1,x,Lpj,Lpr;
1447916e780bShannah_mairs   PetscInt         i,j,nn,r;
1448916e780bShannah_mairs 
1449916e780bShannah_mairs   PetscFunctionBegin;
1450916e780bShannah_mairs   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1451916e780bShannah_mairs   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1452916e780bShannah_mairs   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1453916e780bShannah_mairs 
1454916e780bShannah_mairs   for (j=1; j<p; j++) {
1455916e780bShannah_mairs     x  = gllnodes[j];
1456916e780bShannah_mairs     z0 = 1.;
1457916e780bShannah_mairs     z1 = x;
1458916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1459916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1460916e780bShannah_mairs       z0 = z1;
1461916e780bShannah_mairs       z1 = z2;
1462916e780bShannah_mairs     }
1463916e780bShannah_mairs     Lpj=z2;
1464916e780bShannah_mairs     for (r=1; r<p; r++) {
1465916e780bShannah_mairs       if (r == j) {
1466916e780bShannah_mairs         A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj);
1467916e780bShannah_mairs       } else {
1468916e780bShannah_mairs         x  = gllnodes[r];
1469916e780bShannah_mairs         z0 = 1.;
1470916e780bShannah_mairs         z1 = x;
1471916e780bShannah_mairs         for (nn=1; nn<p; nn++) {
1472916e780bShannah_mairs           z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1473916e780bShannah_mairs           z0 = z1;
1474916e780bShannah_mairs           z1 = z2;
1475916e780bShannah_mairs         }
1476916e780bShannah_mairs         Lpr     = z2;
1477916e780bShannah_mairs         A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r]));
1478916e780bShannah_mairs       }
1479916e780bShannah_mairs     }
1480916e780bShannah_mairs   }
1481916e780bShannah_mairs   for (j=1; j<p+1; j++) {
1482916e780bShannah_mairs     x  = gllnodes[j];
1483916e780bShannah_mairs     z0 = 1.;
1484916e780bShannah_mairs     z1 = x;
1485916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1486916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1487916e780bShannah_mairs       z0 = z1;
1488916e780bShannah_mairs       z1 = z2;
1489916e780bShannah_mairs     }
1490916e780bShannah_mairs     Lpj     = z2;
1491916e780bShannah_mairs     A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j]));
1492916e780bShannah_mairs     A[0][j] = A[j][0];
1493916e780bShannah_mairs   }
1494916e780bShannah_mairs   for (j=0; j<p; j++) {
1495916e780bShannah_mairs     x  = gllnodes[j];
1496916e780bShannah_mairs     z0 = 1.;
1497916e780bShannah_mairs     z1 = x;
1498916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1499916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1500916e780bShannah_mairs       z0 = z1;
1501916e780bShannah_mairs       z1 = z2;
1502916e780bShannah_mairs     }
1503916e780bShannah_mairs     Lpj=z2;
1504916e780bShannah_mairs 
1505916e780bShannah_mairs     A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j]));
1506916e780bShannah_mairs     A[j][p] = A[p][j];
1507916e780bShannah_mairs   }
1508916e780bShannah_mairs   A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.;
1509916e780bShannah_mairs   A[p][p]=A[0][0];
1510916e780bShannah_mairs   *AA = A;
1511916e780bShannah_mairs   PetscFunctionReturn(0);
1512916e780bShannah_mairs }
1513916e780bShannah_mairs 
1514916e780bShannah_mairs /*@C
1515916e780bShannah_mairs    PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element
1516916e780bShannah_mairs 
1517916e780bShannah_mairs    Not Collective
1518916e780bShannah_mairs 
1519916e780bShannah_mairs    Input Parameter:
1520916e780bShannah_mairs +  n - the number of GLL nodes
1521916e780bShannah_mairs .  nodes - the GLL nodes
1522916e780bShannah_mairs .  weights - the GLL weightss
1523916e780bShannah_mairs -  A - the stiffness element
1524916e780bShannah_mairs 
1525916e780bShannah_mairs    Level: beginner
1526916e780bShannah_mairs 
1527916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate()
1528916e780bShannah_mairs 
1529916e780bShannah_mairs @*/
1530916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1531916e780bShannah_mairs {
1532916e780bShannah_mairs   PetscErrorCode ierr;
1533916e780bShannah_mairs 
1534916e780bShannah_mairs   PetscFunctionBegin;
1535916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1536916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1537916e780bShannah_mairs   *AA  = NULL;
1538916e780bShannah_mairs   PetscFunctionReturn(0);
1539916e780bShannah_mairs }
1540916e780bShannah_mairs 
1541916e780bShannah_mairs /*@C
1542916e780bShannah_mairs    PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element
1543916e780bShannah_mairs 
1544916e780bShannah_mairs    Not Collective
1545916e780bShannah_mairs 
1546916e780bShannah_mairs    Input Parameter:
1547916e780bShannah_mairs +  n - the number of GLL nodes
1548916e780bShannah_mairs .  nodes - the GLL nodes
1549916e780bShannah_mairs .  weights - the GLL weights
1550916e780bShannah_mairs 
1551916e780bShannah_mairs    Output Parameter:
1552916e780bShannah_mairs .  AA - the stiffness element
1553916e780bShannah_mairs -  AAT - the transpose of AA (pass in NULL if you do not need this array)
1554916e780bShannah_mairs 
1555916e780bShannah_mairs    Level: beginner
1556916e780bShannah_mairs 
1557916e780bShannah_mairs    Notes:
1558916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementGradientDestroy()
1559916e780bShannah_mairs 
1560916e780bShannah_mairs    You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1561916e780bShannah_mairs 
1562916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1563916e780bShannah_mairs 
1564916e780bShannah_mairs @*/
1565916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
1566916e780bShannah_mairs {
1567916e780bShannah_mairs   PetscReal        **A, **AT = NULL;
1568916e780bShannah_mairs   PetscErrorCode  ierr;
1569916e780bShannah_mairs   const PetscReal  *gllnodes = nodes;
1570916e780bShannah_mairs   const PetscInt   p = n-1;
1571916e780bShannah_mairs   PetscReal        q,qp,Li, Lj,d0;
1572916e780bShannah_mairs   PetscInt         i,j;
1573916e780bShannah_mairs 
1574916e780bShannah_mairs   PetscFunctionBegin;
1575916e780bShannah_mairs   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1576916e780bShannah_mairs   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1577916e780bShannah_mairs   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1578916e780bShannah_mairs 
1579916e780bShannah_mairs   if (AAT) {
1580916e780bShannah_mairs     ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr);
1581916e780bShannah_mairs     ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr);
1582916e780bShannah_mairs     for (i=1; i<n; i++) AT[i] = AT[i-1]+n;
1583916e780bShannah_mairs   }
1584916e780bShannah_mairs 
1585916e780bShannah_mairs   if (n==1) {A[0][0] = 0.;}
1586916e780bShannah_mairs   d0 = (PetscReal)p*((PetscReal)p+1.)/4.;
1587916e780bShannah_mairs   for  (i=0; i<n; i++) {
1588916e780bShannah_mairs     for  (j=0; j<n; j++) {
1589916e780bShannah_mairs       A[i][j] = 0.;
1590fdd31e58Shannah_mairs       qAndLEvaluation(p,gllnodes[i],&q,&qp,&Li);
1591fdd31e58Shannah_mairs       qAndLEvaluation(p,gllnodes[j],&q,&qp,&Lj);
1592916e780bShannah_mairs       if (i!=j)             A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j]));
1593916e780bShannah_mairs       if ((j==i) && (i==0)) A[i][j] = -d0;
1594916e780bShannah_mairs       if (j==i && i==p)     A[i][j] = d0;
1595916e780bShannah_mairs       if (AT) AT[j][i] = A[i][j];
1596916e780bShannah_mairs     }
1597916e780bShannah_mairs   }
1598916e780bShannah_mairs   if (AAT) *AAT = AT;
1599916e780bShannah_mairs   *AA  = A;
1600916e780bShannah_mairs   PetscFunctionReturn(0);
1601916e780bShannah_mairs }
1602916e780bShannah_mairs 
1603916e780bShannah_mairs /*@C
1604916e780bShannah_mairs    PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate()
1605916e780bShannah_mairs 
1606916e780bShannah_mairs    Not Collective
1607916e780bShannah_mairs 
1608916e780bShannah_mairs    Input Parameter:
1609916e780bShannah_mairs +  n - the number of GLL nodes
1610916e780bShannah_mairs .  nodes - the GLL nodes
1611916e780bShannah_mairs .  weights - the GLL weights
1612916e780bShannah_mairs .  AA - the stiffness element
1613916e780bShannah_mairs -  AAT - the transpose of the element
1614916e780bShannah_mairs 
1615916e780bShannah_mairs    Level: beginner
1616916e780bShannah_mairs 
1617916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate()
1618916e780bShannah_mairs 
1619916e780bShannah_mairs @*/
1620916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
1621916e780bShannah_mairs {
1622916e780bShannah_mairs   PetscErrorCode ierr;
1623916e780bShannah_mairs 
1624916e780bShannah_mairs   PetscFunctionBegin;
1625916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1626916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1627916e780bShannah_mairs   *AA  = NULL;
1628916e780bShannah_mairs   if (*AAT) {
1629916e780bShannah_mairs     ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr);
1630916e780bShannah_mairs     ierr = PetscFree(*AAT);CHKERRQ(ierr);
1631916e780bShannah_mairs     *AAT  = NULL;
1632916e780bShannah_mairs   }
1633916e780bShannah_mairs   PetscFunctionReturn(0);
1634916e780bShannah_mairs }
1635916e780bShannah_mairs 
1636916e780bShannah_mairs /*@C
1637916e780bShannah_mairs    PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element
1638916e780bShannah_mairs 
1639916e780bShannah_mairs    Not Collective
1640916e780bShannah_mairs 
1641916e780bShannah_mairs    Input Parameter:
1642916e780bShannah_mairs +  n - the number of GLL nodes
1643916e780bShannah_mairs .  nodes - the GLL nodes
1644916e780bShannah_mairs .  weights - the GLL weightss
1645916e780bShannah_mairs 
1646916e780bShannah_mairs    Output Parameter:
1647916e780bShannah_mairs .  AA - the stiffness element
1648916e780bShannah_mairs 
1649916e780bShannah_mairs    Level: beginner
1650916e780bShannah_mairs 
1651916e780bShannah_mairs    Notes:
1652916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy()
1653916e780bShannah_mairs 
1654916e780bShannah_mairs    This is the same as the Gradient operator multiplied by the diagonal mass matrix
1655916e780bShannah_mairs 
1656916e780bShannah_mairs    You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1657916e780bShannah_mairs 
1658916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy()
1659916e780bShannah_mairs 
1660916e780bShannah_mairs @*/
1661916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1662916e780bShannah_mairs {
1663916e780bShannah_mairs   PetscReal       **D;
1664916e780bShannah_mairs   PetscErrorCode  ierr;
1665916e780bShannah_mairs   const PetscReal  *gllweights = weights;
1666916e780bShannah_mairs   const PetscInt   glln = n;
1667916e780bShannah_mairs   PetscInt         i,j;
1668916e780bShannah_mairs 
1669916e780bShannah_mairs   PetscFunctionBegin;
1670916e780bShannah_mairs   ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr);
1671916e780bShannah_mairs   for (i=0; i<glln; i++){
1672916e780bShannah_mairs     for (j=0; j<glln; j++) {
1673916e780bShannah_mairs       D[i][j] = gllweights[i]*D[i][j];
1674916e780bShannah_mairs     }
1675916e780bShannah_mairs   }
1676916e780bShannah_mairs   *AA = D;
1677916e780bShannah_mairs   PetscFunctionReturn(0);
1678916e780bShannah_mairs }
1679916e780bShannah_mairs 
1680916e780bShannah_mairs /*@C
1681916e780bShannah_mairs    PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element
1682916e780bShannah_mairs 
1683916e780bShannah_mairs    Not Collective
1684916e780bShannah_mairs 
1685916e780bShannah_mairs    Input Parameter:
1686916e780bShannah_mairs +  n - the number of GLL nodes
1687916e780bShannah_mairs .  nodes - the GLL nodes
1688916e780bShannah_mairs .  weights - the GLL weights
1689916e780bShannah_mairs -  A - advection
1690916e780bShannah_mairs 
1691916e780bShannah_mairs    Level: beginner
1692916e780bShannah_mairs 
1693916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate()
1694916e780bShannah_mairs 
1695916e780bShannah_mairs @*/
1696916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1697916e780bShannah_mairs {
1698916e780bShannah_mairs   PetscErrorCode ierr;
1699916e780bShannah_mairs 
1700916e780bShannah_mairs   PetscFunctionBegin;
1701916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1702916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1703916e780bShannah_mairs   *AA  = NULL;
1704916e780bShannah_mairs   PetscFunctionReturn(0);
1705916e780bShannah_mairs }
1706916e780bShannah_mairs 
1707916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1708916e780bShannah_mairs {
1709916e780bShannah_mairs   PetscReal        **A;
1710916e780bShannah_mairs   PetscErrorCode  ierr;
1711916e780bShannah_mairs   const PetscReal  *gllweights = weights;
1712916e780bShannah_mairs   const PetscInt   glln = n;
1713916e780bShannah_mairs   PetscInt         i,j;
1714916e780bShannah_mairs 
1715916e780bShannah_mairs   PetscFunctionBegin;
1716916e780bShannah_mairs   ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr);
1717916e780bShannah_mairs   ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr);
1718916e780bShannah_mairs   for (i=1; i<glln; i++) A[i] = A[i-1]+glln;
1719916e780bShannah_mairs   if (glln==1) {A[0][0] = 0.;}
1720916e780bShannah_mairs   for  (i=0; i<glln; i++) {
1721916e780bShannah_mairs     for  (j=0; j<glln; j++) {
1722916e780bShannah_mairs       A[i][j] = 0.;
1723916e780bShannah_mairs       if (j==i)     A[i][j] = gllweights[i];
1724916e780bShannah_mairs     }
1725916e780bShannah_mairs   }
1726916e780bShannah_mairs   *AA  = A;
1727916e780bShannah_mairs   PetscFunctionReturn(0);
1728916e780bShannah_mairs }
1729916e780bShannah_mairs 
1730916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1731916e780bShannah_mairs {
1732916e780bShannah_mairs   PetscErrorCode ierr;
1733916e780bShannah_mairs 
1734916e780bShannah_mairs   PetscFunctionBegin;
1735916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1736916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1737916e780bShannah_mairs   *AA  = NULL;
1738916e780bShannah_mairs   PetscFunctionReturn(0);
1739916e780bShannah_mairs }
1740916e780bShannah_mairs 
1741