xref: /petsc/src/dm/dt/interface/dt.c (revision 916e780b13b6ed67b2ef5efd8e6dd7dbdbb7e8c1)
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;
314d9bac1caSLisandro Dalcin   if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "Quadrature on %D points with %D components of order %D\n", quad->numPoints, quad->Nc, quad->order);CHKERRQ(ierr);}
315d9bac1caSLisandro Dalcin   else              {ierr = PetscViewerASCIIPrintf(v, "Quadrature on %D points of order %D\n", quad->numPoints, quad->order);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) {
319d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(v, "(");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);
326d9bac1caSLisandro Dalcin     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "(");CHKERRQ(ierr);}
327d9bac1caSLisandro Dalcin     for (c = 0; c < quad->Nc; ++c) {
328d9bac1caSLisandro Dalcin       if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
329d9bac1caSLisandro Dalcin       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);}
36198c3331eSBarry Smith   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad, viewer);CHKERRQ(ierr);
362d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
363d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
364d9bac1caSLisandro Dalcin   if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);}
365d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
366bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
367bfa639d9SMatthew G. Knepley }
368bfa639d9SMatthew G. Knepley 
36989710940SMatthew G. Knepley /*@C
37089710940SMatthew G. Knepley   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
37189710940SMatthew G. Knepley 
37289710940SMatthew G. Knepley   Not collective
37389710940SMatthew G. Knepley 
37489710940SMatthew G. Knepley   Input Parameter:
37589710940SMatthew G. Knepley + q - The original PetscQuadrature
37689710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into
37789710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement
37889710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement
37989710940SMatthew G. Knepley 
38089710940SMatthew G. Knepley   Output Parameters:
38189710940SMatthew G. Knepley . dim - The dimension
38289710940SMatthew G. Knepley 
38389710940SMatthew G. Knepley   Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
38489710940SMatthew G. Knepley 
385f5f57ec0SBarry Smith  Not available from Fortran
386f5f57ec0SBarry Smith 
38789710940SMatthew G. Knepley   Level: intermediate
38889710940SMatthew G. Knepley 
38989710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
39089710940SMatthew G. Knepley @*/
39189710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
39289710940SMatthew G. Knepley {
39389710940SMatthew G. Knepley   const PetscReal *points,    *weights;
39489710940SMatthew G. Knepley   PetscReal       *pointsRef, *weightsRef;
395a6b92713SMatthew G. Knepley   PetscInt         dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
39689710940SMatthew G. Knepley   PetscErrorCode   ierr;
39789710940SMatthew G. Knepley 
39889710940SMatthew G. Knepley   PetscFunctionBegin;
39989710940SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
40089710940SMatthew G. Knepley   PetscValidPointer(v0, 3);
40189710940SMatthew G. Knepley   PetscValidPointer(jac, 4);
40289710940SMatthew G. Knepley   PetscValidPointer(qref, 5);
40389710940SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr);
40489710940SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
405a6b92713SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr);
40689710940SMatthew G. Knepley   npointsRef = npoints*numSubelements;
40789710940SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr);
408a6b92713SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr);
40989710940SMatthew G. Knepley   for (c = 0; c < numSubelements; ++c) {
41089710940SMatthew G. Knepley     for (p = 0; p < npoints; ++p) {
41189710940SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
41289710940SMatthew G. Knepley         pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
41389710940SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
41489710940SMatthew G. Knepley           pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
41589710940SMatthew G. Knepley         }
41689710940SMatthew G. Knepley       }
41789710940SMatthew G. Knepley       /* Could also use detJ here */
418a6b92713SMatthew G. Knepley       for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements;
41989710940SMatthew G. Knepley     }
42089710940SMatthew G. Knepley   }
42189710940SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr);
422a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr);
42389710940SMatthew G. Knepley   PetscFunctionReturn(0);
42489710940SMatthew G. Knepley }
42589710940SMatthew G. Knepley 
42637045ce4SJed Brown /*@
42737045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
42837045ce4SJed Brown 
42937045ce4SJed Brown    Not Collective
43037045ce4SJed Brown 
43137045ce4SJed Brown    Input Arguments:
43237045ce4SJed Brown +  npoints - number of spatial points to evaluate at
43337045ce4SJed Brown .  points - array of locations to evaluate at
43437045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
43537045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
43637045ce4SJed Brown 
43737045ce4SJed Brown    Output Arguments:
4380298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
4390298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
4400298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
44137045ce4SJed Brown 
44237045ce4SJed Brown    Level: intermediate
44337045ce4SJed Brown 
44437045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
44537045ce4SJed Brown @*/
44637045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
44737045ce4SJed Brown {
44837045ce4SJed Brown   PetscInt i,maxdegree;
44937045ce4SJed Brown 
45037045ce4SJed Brown   PetscFunctionBegin;
45137045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
45237045ce4SJed Brown   maxdegree = degrees[ndegree-1];
45337045ce4SJed Brown   for (i=0; i<npoints; i++) {
45437045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
45537045ce4SJed Brown     PetscInt  j,k;
45637045ce4SJed Brown     x    = points[i];
45737045ce4SJed Brown     pm2  = 0;
45837045ce4SJed Brown     pm1  = 1;
45937045ce4SJed Brown     pd2  = 0;
46037045ce4SJed Brown     pd1  = 0;
46137045ce4SJed Brown     pdd2 = 0;
46237045ce4SJed Brown     pdd1 = 0;
46337045ce4SJed Brown     k    = 0;
46437045ce4SJed Brown     if (degrees[k] == 0) {
46537045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
46637045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
46737045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
46837045ce4SJed Brown       k++;
46937045ce4SJed Brown     }
47037045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
47137045ce4SJed Brown       PetscReal p,d,dd;
47237045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
47337045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
47437045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
47537045ce4SJed Brown       pm2  = pm1;
47637045ce4SJed Brown       pm1  = p;
47737045ce4SJed Brown       pd2  = pd1;
47837045ce4SJed Brown       pd1  = d;
47937045ce4SJed Brown       pdd2 = pdd1;
48037045ce4SJed Brown       pdd1 = dd;
48137045ce4SJed Brown       if (degrees[k] == j) {
48237045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
48337045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
48437045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
48537045ce4SJed Brown       }
48637045ce4SJed Brown     }
48737045ce4SJed Brown   }
48837045ce4SJed Brown   PetscFunctionReturn(0);
48937045ce4SJed Brown }
49037045ce4SJed Brown 
49137045ce4SJed Brown /*@
49237045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
49337045ce4SJed Brown 
49437045ce4SJed Brown    Not Collective
49537045ce4SJed Brown 
49637045ce4SJed Brown    Input Arguments:
49737045ce4SJed Brown +  npoints - number of points
49837045ce4SJed Brown .  a - left end of interval (often-1)
49937045ce4SJed Brown -  b - right end of interval (often +1)
50037045ce4SJed Brown 
50137045ce4SJed Brown    Output Arguments:
50237045ce4SJed Brown +  x - quadrature points
50337045ce4SJed Brown -  w - quadrature weights
50437045ce4SJed Brown 
50537045ce4SJed Brown    Level: intermediate
50637045ce4SJed Brown 
50737045ce4SJed Brown    References:
50896a0c994SBarry Smith .   1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
50937045ce4SJed Brown 
51037045ce4SJed Brown .seealso: PetscDTLegendreEval()
51137045ce4SJed Brown @*/
51237045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
51337045ce4SJed Brown {
51437045ce4SJed Brown   PetscErrorCode ierr;
51537045ce4SJed Brown   PetscInt       i;
51637045ce4SJed Brown   PetscReal      *work;
51737045ce4SJed Brown   PetscScalar    *Z;
51837045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
51937045ce4SJed Brown 
52037045ce4SJed Brown   PetscFunctionBegin;
5210bfcf5a5SMatthew G. Knepley   ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr);
52237045ce4SJed Brown   /* Set up the Golub-Welsch system */
52337045ce4SJed Brown   for (i=0; i<npoints; i++) {
52437045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
52537045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
52637045ce4SJed Brown   }
527dcca6d9dSJed Brown   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
528c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
52937045ce4SJed Brown   LDZ  = N;
53037045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
5318b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
53237045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
5331c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
53437045ce4SJed Brown 
53537045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
53637045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
53737045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
53819a57d60SBarry Smith     if (x[i] == -0.0) x[i] = 0.0;
53937045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
5400d644c17SKarl Rupp 
54188393a60SJed Brown     w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
54237045ce4SJed Brown   }
54337045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
54437045ce4SJed Brown   PetscFunctionReturn(0);
54537045ce4SJed Brown }
546194825f6SJed Brown 
5478272889dSSatish Balay static void qAndLEvaluation(PetscInt n, PetscReal x, PetscReal *q, PetscReal *qp, PetscReal *Ln)
5488272889dSSatish Balay /*
5498272889dSSatish Balay   Compute the polynomial q(x) = L_{N+1}(x) - L_{n-1}(x) and its derivative in
5508272889dSSatish Balay   addition to L_N(x) as these are needed for computing the GLL points via Newton's method.
5518272889dSSatish Balay   Reference: "Implementing Spectral Methods for Partial Differential Equations: Algorithms
5528272889dSSatish Balay   for Scientists and Engineers" by David A. Kopriva.
5538272889dSSatish Balay */
5548272889dSSatish Balay {
5558272889dSSatish Balay   PetscInt k;
5568272889dSSatish Balay 
5578272889dSSatish Balay   PetscReal Lnp;
5588272889dSSatish Balay   PetscReal Lnp1, Lnp1p;
5598272889dSSatish Balay   PetscReal Lnm1, Lnm1p;
5608272889dSSatish Balay   PetscReal Lnm2, Lnm2p;
5618272889dSSatish Balay 
5628272889dSSatish Balay   Lnm1  = 1.0;
5638272889dSSatish Balay   *Ln   = x;
5648272889dSSatish Balay   Lnm1p = 0.0;
5658272889dSSatish Balay   Lnp   = 1.0;
5668272889dSSatish Balay 
5678272889dSSatish Balay   for (k=2; k<=n; ++k) {
5688272889dSSatish Balay     Lnm2  = Lnm1;
5698272889dSSatish Balay     Lnm1  = *Ln;
5708272889dSSatish Balay     Lnm2p = Lnm1p;
5718272889dSSatish Balay     Lnm1p = Lnp;
5728272889dSSatish Balay     *Ln   = (2.*((PetscReal)k)-1.)/(1.0*((PetscReal)k))*x*Lnm1 - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm2;
5738272889dSSatish Balay     Lnp   = Lnm2p + (2.0*((PetscReal)k)-1.)*Lnm1;
5748272889dSSatish Balay   }
5758272889dSSatish Balay   k     = n+1;
5768272889dSSatish Balay   Lnp1  = (2.*((PetscReal)k)-1.)/(((PetscReal)k))*x*(*Ln) - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm1;
5778272889dSSatish Balay   Lnp1p = Lnm1p + (2.0*((PetscReal)k)-1.)*(*Ln);
5788272889dSSatish Balay   *q    = Lnp1 - Lnm1;
5798272889dSSatish Balay   *qp   = Lnp1p - Lnm1p;
5808272889dSSatish Balay }
5818272889dSSatish Balay 
5828272889dSSatish Balay /*@C
5838272889dSSatish Balay    PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre
5848272889dSSatish Balay                       nodes of a given size on the domain [-1,1]
5858272889dSSatish Balay 
5868272889dSSatish Balay    Not Collective
5878272889dSSatish Balay 
5888272889dSSatish Balay    Input Parameter:
5898272889dSSatish Balay +  n - number of grid nodes
590*916e780bShannah_mairs -  type - PETSCGaussLobattoLegendre_VIA_LINEARALGEBRA or PETSCGaussLobattoLegendre_VIA_NEWTON
5918272889dSSatish Balay 
5928272889dSSatish Balay    Output Arguments:
5938272889dSSatish Balay +  x - quadrature points
5948272889dSSatish Balay -  w - quadrature weights
5958272889dSSatish Balay 
5968272889dSSatish Balay    Notes:
5978272889dSSatish Balay     For n > 30  the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
5988272889dSSatish Balay           close enough to the desired solution
5998272889dSSatish Balay 
6008272889dSSatish Balay    These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes
6018272889dSSatish Balay 
6028272889dSSatish 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
6038272889dSSatish Balay 
6048272889dSSatish Balay    Level: intermediate
6058272889dSSatish Balay 
6068272889dSSatish Balay .seealso: PetscDTGaussQuadrature()
6078272889dSSatish Balay 
6088272889dSSatish Balay @*/
609*916e780bShannah_mairs PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w)
6108272889dSSatish Balay {
6118272889dSSatish Balay   PetscErrorCode ierr;
6128272889dSSatish Balay 
6138272889dSSatish Balay   PetscFunctionBegin;
6148272889dSSatish Balay   if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element");
6158272889dSSatish Balay 
616*916e780bShannah_mairs   if (type == PETSCGaussLobattoLegendre_VIA_LINEARALGEBRA) {
6178272889dSSatish Balay     PetscReal      *M,si;
6188272889dSSatish Balay     PetscBLASInt   bn,lierr;
6198272889dSSatish Balay     PetscReal      x0,z0,z1,z2;
6208272889dSSatish Balay     PetscInt       i,p = npoints - 1,nn;
6218272889dSSatish Balay 
6228272889dSSatish Balay     x[0]   =-1.0;
6238272889dSSatish Balay     x[npoints-1] = 1.0;
6248272889dSSatish Balay     if (npoints-2 > 0){
6258272889dSSatish Balay       ierr = PetscMalloc1(npoints-1,&M);CHKERRQ(ierr);
6268272889dSSatish Balay       for (i=0; i<npoints-2; i++) {
6278272889dSSatish Balay         si  = ((PetscReal)i)+1.0;
6288272889dSSatish Balay         M[i]=0.5*PetscSqrtReal(si*(si+2.0)/((si+0.5)*(si+1.5)));
6298272889dSSatish Balay       }
6308272889dSSatish Balay       ierr = PetscBLASIntCast(npoints-2,&bn);CHKERRQ(ierr);
6318272889dSSatish Balay       ierr = PetscMemzero(&x[1],bn*sizeof(x[1]));CHKERRQ(ierr);
6328272889dSSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
6338272889dSSatish Balay       x0=0;
6348272889dSSatish Balay       PetscStackCallBLAS("LAPACKsteqr",LAPACKREALsteqr_("N",&bn,&x[1],M,&x0,&bn,M,&lierr));
6358272889dSSatish Balay       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in STERF Lapack routine %d",(int)lierr);
6368272889dSSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
6378272889dSSatish Balay       ierr = PetscFree(M);CHKERRQ(ierr);
6388272889dSSatish Balay     }
6398272889dSSatish Balay     if ((npoints-1)%2==0) {
6408272889dSSatish Balay       x[(npoints-1)/2]   = 0.0; /* hard wire to exactly 0.0 since linear algebra produces nonzero */
6418272889dSSatish Balay     }
6428272889dSSatish Balay 
6438272889dSSatish Balay     w[0] = w[p] = 2.0/(((PetscReal)(p))*(((PetscReal)p)+1.0));
6448272889dSSatish Balay     z2 = -1.;                      /* Dummy value to avoid -Wmaybe-initialized */
6458272889dSSatish Balay     for (i=1; i<p; i++) {
6468272889dSSatish Balay       x0  = x[i];
6478272889dSSatish Balay       z0 = 1.0;
6488272889dSSatish Balay       z1 = x0;
6498272889dSSatish Balay       for (nn=1; nn<p; nn++) {
6508272889dSSatish Balay         z2 = x0*z1*(2.0*((PetscReal)nn)+1.0)/(((PetscReal)nn)+1.0)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.0));
6518272889dSSatish Balay         z0 = z1;
6528272889dSSatish Balay         z1 = z2;
6538272889dSSatish Balay       }
6548272889dSSatish Balay       w[i]=2.0/(((PetscReal)p)*(((PetscReal)p)+1.0)*z2*z2);
6558272889dSSatish Balay     }
6568272889dSSatish Balay   } else {
6578272889dSSatish Balay     PetscInt  j,m;
6588272889dSSatish Balay     PetscReal z1,z,q,qp,Ln;
6598272889dSSatish Balay     PetscReal *pt;
6608272889dSSatish Balay     ierr = PetscMalloc1(npoints,&pt);CHKERRQ(ierr);
6618272889dSSatish Balay 
662*916e780bShannah_mairs     if (npoints > 30) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PETSCGaussLobattoLegendre_VIA_NEWTON produces incorrect answers for n > 30");
6638272889dSSatish Balay     x[0]     = -1.0;
6648272889dSSatish Balay     x[npoints-1]   = 1.0;
6658272889dSSatish Balay     w[0]   = w[npoints-1] = 2./(((PetscReal)npoints)*(((PetscReal)npoints)-1.0));;
6668272889dSSatish Balay     m  = (npoints-1)/2; /* The roots are symmetric, so we only find half of them. */
6678272889dSSatish Balay     for (j=1; j<=m; j++) { /* Loop over the desired roots. */
6688272889dSSatish 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));
6698272889dSSatish Balay       /* Starting with the above approximation to the ith root, we enter */
6708272889dSSatish Balay       /* the main loop of refinement by Newton's method.                 */
6718272889dSSatish Balay       do {
6728272889dSSatish Balay         qAndLEvaluation(npoints-1,z,&q,&qp,&Ln);
6738272889dSSatish Balay         z1 = z;
6748272889dSSatish Balay         z  = z1-q/qp; /* Newton's method. */
6758272889dSSatish Balay       } while (PetscAbs(z-z1) > 10.*PETSC_MACHINE_EPSILON);
6768272889dSSatish Balay       qAndLEvaluation(npoints-1,z,&q,&qp,&Ln);
6778272889dSSatish Balay 
6788272889dSSatish Balay       x[j]       = z;
6798272889dSSatish Balay       x[npoints-1-j]   = -z;      /* and put in its symmetric counterpart.   */
6808272889dSSatish Balay       w[j]     = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln);  /* Compute the weight */
6818272889dSSatish Balay       w[npoints-1-j] = w[j];                 /* and its symmetric counterpart. */
6828272889dSSatish Balay       pt[j]=qp;
6838272889dSSatish Balay     }
6848272889dSSatish Balay 
6858272889dSSatish Balay     if ((npoints-1)%2==0) {
6868272889dSSatish Balay       qAndLEvaluation(npoints-1,0.0,&q,&qp,&Ln);
6878272889dSSatish Balay       x[(npoints-1)/2]   = 0.0;
6888272889dSSatish Balay       w[(npoints-1)/2] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln);
6898272889dSSatish Balay     }
6908272889dSSatish Balay     ierr = PetscFree(pt);CHKERRQ(ierr);
6918272889dSSatish Balay   }
6928272889dSSatish Balay   PetscFunctionReturn(0);
6938272889dSSatish Balay }
6948272889dSSatish Balay 
695744bafbcSMatthew G. Knepley /*@
696744bafbcSMatthew G. Knepley   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
697744bafbcSMatthew G. Knepley 
698744bafbcSMatthew G. Knepley   Not Collective
699744bafbcSMatthew G. Knepley 
700744bafbcSMatthew G. Knepley   Input Arguments:
701744bafbcSMatthew G. Knepley + dim     - The spatial dimension
702a6b92713SMatthew G. Knepley . Nc      - The number of components
703744bafbcSMatthew G. Knepley . npoints - number of points in one dimension
704744bafbcSMatthew G. Knepley . a       - left end of interval (often-1)
705744bafbcSMatthew G. Knepley - b       - right end of interval (often +1)
706744bafbcSMatthew G. Knepley 
707744bafbcSMatthew G. Knepley   Output Argument:
708744bafbcSMatthew G. Knepley . q - A PetscQuadrature object
709744bafbcSMatthew G. Knepley 
710744bafbcSMatthew G. Knepley   Level: intermediate
711744bafbcSMatthew G. Knepley 
712744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
713744bafbcSMatthew G. Knepley @*/
714a6b92713SMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
715744bafbcSMatthew G. Knepley {
716a6b92713SMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
717744bafbcSMatthew G. Knepley   PetscReal     *x, *w, *xw, *ww;
718744bafbcSMatthew G. Knepley   PetscErrorCode ierr;
719744bafbcSMatthew G. Knepley 
720744bafbcSMatthew G. Knepley   PetscFunctionBegin;
721744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
722a6b92713SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr);
723744bafbcSMatthew G. Knepley   /* Set up the Golub-Welsch system */
724744bafbcSMatthew G. Knepley   switch (dim) {
725744bafbcSMatthew G. Knepley   case 0:
726744bafbcSMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
727744bafbcSMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
728744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
729a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
730744bafbcSMatthew G. Knepley     x[0] = 0.0;
731a6b92713SMatthew G. Knepley     for (c = 0; c < Nc; ++c) w[c] = 1.0;
732744bafbcSMatthew G. Knepley     break;
733744bafbcSMatthew G. Knepley   case 1:
734a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr);
735a6b92713SMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr);
736a6b92713SMatthew G. Knepley     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i];
737a6b92713SMatthew G. Knepley     ierr = PetscFree(ww);CHKERRQ(ierr);
738744bafbcSMatthew G. Knepley     break;
739744bafbcSMatthew G. Knepley   case 2:
740744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
741744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
742744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
743744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
744744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+0] = xw[i];
745744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+1] = xw[j];
746a6b92713SMatthew G. Knepley         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j];
747744bafbcSMatthew G. Knepley       }
748744bafbcSMatthew G. Knepley     }
749744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
750744bafbcSMatthew G. Knepley     break;
751744bafbcSMatthew G. Knepley   case 3:
752744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
753744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
754744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
755744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
756744bafbcSMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
757744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
758744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
759744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
760a6b92713SMatthew G. Knepley           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k];
761744bafbcSMatthew G. Knepley         }
762744bafbcSMatthew G. Knepley       }
763744bafbcSMatthew G. Knepley     }
764744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
765744bafbcSMatthew G. Knepley     break;
766744bafbcSMatthew G. Knepley   default:
767744bafbcSMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
768744bafbcSMatthew G. Knepley   }
769744bafbcSMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
7702f5fb066SToby Isaac   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
771a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
772d9bac1caSLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr);
773744bafbcSMatthew G. Knepley   PetscFunctionReturn(0);
774744bafbcSMatthew G. Knepley }
775744bafbcSMatthew G. Knepley 
776494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
777494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
778494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
779494e7359SMatthew G. Knepley {
780494e7359SMatthew G. Knepley   PetscReal f = 1.0;
781494e7359SMatthew G. Knepley   PetscInt  i;
782494e7359SMatthew G. Knepley 
783494e7359SMatthew G. Knepley   PetscFunctionBegin;
784494e7359SMatthew G. Knepley   for (i = 1; i < n+1; ++i) f *= i;
785494e7359SMatthew G. Knepley   *factorial = f;
786494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
787494e7359SMatthew G. Knepley }
788494e7359SMatthew G. Knepley 
789494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
790494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
791494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
792494e7359SMatthew G. Knepley {
793494e7359SMatthew G. Knepley   PetscReal apb, pn1, pn2;
794494e7359SMatthew G. Knepley   PetscInt  k;
795494e7359SMatthew G. Knepley 
796494e7359SMatthew G. Knepley   PetscFunctionBegin;
797494e7359SMatthew G. Knepley   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
798494e7359SMatthew G. Knepley   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
799494e7359SMatthew G. Knepley   apb = a + b;
800494e7359SMatthew G. Knepley   pn2 = 1.0;
801494e7359SMatthew G. Knepley   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
802494e7359SMatthew G. Knepley   *P  = 0.0;
803494e7359SMatthew G. Knepley   for (k = 2; k < n+1; ++k) {
804494e7359SMatthew G. Knepley     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
805494e7359SMatthew G. Knepley     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
806494e7359SMatthew G. Knepley     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
807494e7359SMatthew G. Knepley     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
808494e7359SMatthew G. Knepley 
809494e7359SMatthew G. Knepley     a2  = a2 / a1;
810494e7359SMatthew G. Knepley     a3  = a3 / a1;
811494e7359SMatthew G. Knepley     a4  = a4 / a1;
812494e7359SMatthew G. Knepley     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
813494e7359SMatthew G. Knepley     pn2 = pn1;
814494e7359SMatthew G. Knepley     pn1 = *P;
815494e7359SMatthew G. Knepley   }
816494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
817494e7359SMatthew G. Knepley }
818494e7359SMatthew G. Knepley 
819494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
820494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
821494e7359SMatthew G. Knepley {
822494e7359SMatthew G. Knepley   PetscReal      nP;
823494e7359SMatthew G. Knepley   PetscErrorCode ierr;
824494e7359SMatthew G. Knepley 
825494e7359SMatthew G. Knepley   PetscFunctionBegin;
826494e7359SMatthew G. Knepley   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
827494e7359SMatthew G. Knepley   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
828494e7359SMatthew G. Knepley   *P   = 0.5 * (a + b + n + 1) * nP;
829494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
830494e7359SMatthew G. Knepley }
831494e7359SMatthew G. Knepley 
832494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
833494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
834494e7359SMatthew G. Knepley {
835494e7359SMatthew G. Knepley   PetscFunctionBegin;
836494e7359SMatthew G. Knepley   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
837494e7359SMatthew G. Knepley   *eta = y;
838494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
839494e7359SMatthew G. Knepley }
840494e7359SMatthew G. Knepley 
841494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
842494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
843494e7359SMatthew G. Knepley {
844494e7359SMatthew G. Knepley   PetscFunctionBegin;
845494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
846494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
847494e7359SMatthew G. Knepley   *zeta = z;
848494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
849494e7359SMatthew G. Knepley }
850494e7359SMatthew G. Knepley 
851494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
852494e7359SMatthew G. Knepley {
853494e7359SMatthew G. Knepley   PetscInt       maxIter = 100;
854494e7359SMatthew G. Knepley   PetscReal      eps     = 1.0e-8;
855a8291ba1SSatish Balay   PetscReal      a1, a2, a3, a4, a5, a6;
856494e7359SMatthew G. Knepley   PetscInt       k;
857494e7359SMatthew G. Knepley   PetscErrorCode ierr;
858494e7359SMatthew G. Knepley 
859494e7359SMatthew G. Knepley   PetscFunctionBegin;
860a8291ba1SSatish Balay 
8618b49ba18SBarry Smith   a1      = PetscPowReal(2.0, a+b+1);
862a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA)
8630646a658SBarry Smith   a2      = PetscTGamma(a + npoints + 1);
8640646a658SBarry Smith   a3      = PetscTGamma(b + npoints + 1);
8650646a658SBarry Smith   a4      = PetscTGamma(a + b + npoints + 1);
866a8291ba1SSatish Balay #else
86729bcbfd0SToby Isaac   {
868d24bbb91SToby Isaac     PetscInt ia, ib;
86929bcbfd0SToby Isaac 
870d24bbb91SToby Isaac     ia = (PetscInt) a;
871d24bbb91SToby Isaac     ib = (PetscInt) b;
872d24bbb91SToby 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 */
873d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr);
874d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr);
875d24bbb91SToby Isaac       ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr);
87629bcbfd0SToby Isaac     } else {
877a8291ba1SSatish Balay       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
87829bcbfd0SToby Isaac     }
87929bcbfd0SToby Isaac   }
880a8291ba1SSatish Balay #endif
881a8291ba1SSatish Balay 
882494e7359SMatthew G. Knepley   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
883494e7359SMatthew G. Knepley   a6   = a1 * a2 * a3 / a4 / a5;
884494e7359SMatthew G. Knepley   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
885494e7359SMatthew G. Knepley    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
886494e7359SMatthew G. Knepley   for (k = 0; k < npoints; ++k) {
8878b49ba18SBarry Smith     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
888494e7359SMatthew G. Knepley     PetscInt  j;
889494e7359SMatthew G. Knepley 
890494e7359SMatthew G. Knepley     if (k > 0) r = 0.5 * (r + x[k-1]);
891494e7359SMatthew G. Knepley     for (j = 0; j < maxIter; ++j) {
892494e7359SMatthew G. Knepley       PetscReal s = 0.0, delta, f, fp;
893494e7359SMatthew G. Knepley       PetscInt  i;
894494e7359SMatthew G. Knepley 
895494e7359SMatthew G. Knepley       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
896494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
897494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
898494e7359SMatthew G. Knepley       delta = f / (fp - f * s);
899494e7359SMatthew G. Knepley       r     = r - delta;
90077b4d14cSPeter Brune       if (PetscAbsReal(delta) < eps) break;
901494e7359SMatthew G. Knepley     }
902494e7359SMatthew G. Knepley     x[k] = r;
903494e7359SMatthew G. Knepley     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
904494e7359SMatthew G. Knepley     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
905494e7359SMatthew G. Knepley   }
906494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
907494e7359SMatthew G. Knepley }
908494e7359SMatthew G. Knepley 
909f5f57ec0SBarry Smith /*@
910494e7359SMatthew G. Knepley   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
911494e7359SMatthew G. Knepley 
912494e7359SMatthew G. Knepley   Not Collective
913494e7359SMatthew G. Knepley 
914494e7359SMatthew G. Knepley   Input Arguments:
915494e7359SMatthew G. Knepley + dim     - The simplex dimension
916a6b92713SMatthew G. Knepley . Nc      - The number of components
917dcce0ee2SMatthew G. Knepley . npoints - The number of points in one dimension
918494e7359SMatthew G. Knepley . a       - left end of interval (often-1)
919494e7359SMatthew G. Knepley - b       - right end of interval (often +1)
920494e7359SMatthew G. Knepley 
921744bafbcSMatthew G. Knepley   Output Argument:
922552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
923494e7359SMatthew G. Knepley 
924494e7359SMatthew G. Knepley   Level: intermediate
925494e7359SMatthew G. Knepley 
926494e7359SMatthew G. Knepley   References:
92796a0c994SBarry Smith .  1. - Karniadakis and Sherwin.  FIAT
928494e7359SMatthew G. Knepley 
929744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
930494e7359SMatthew G. Knepley @*/
931dcce0ee2SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
932494e7359SMatthew G. Knepley {
933dcce0ee2SMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints;
934494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
935a6b92713SMatthew G. Knepley   PetscInt       i, j, k, c;
936494e7359SMatthew G. Knepley   PetscErrorCode ierr;
937494e7359SMatthew G. Knepley 
938494e7359SMatthew G. Knepley   PetscFunctionBegin;
939494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
940dcce0ee2SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr);
941dcce0ee2SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr);
942494e7359SMatthew G. Knepley   switch (dim) {
943707aa5c5SMatthew G. Knepley   case 0:
944707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
945707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
946785e854fSJed Brown     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
947a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
948707aa5c5SMatthew G. Knepley     x[0] = 0.0;
949a6b92713SMatthew G. Knepley     for (c = 0; c < Nc; ++c) w[c] = 1.0;
950707aa5c5SMatthew G. Knepley     break;
951494e7359SMatthew G. Knepley   case 1:
952dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr);
953dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr);
954dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i];
955a6b92713SMatthew G. Knepley     ierr = PetscFree(wx);CHKERRQ(ierr);
956494e7359SMatthew G. Knepley     break;
957494e7359SMatthew G. Knepley   case 2:
958dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr);
959dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
960dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
961dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
962dcce0ee2SMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
963dcce0ee2SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr);
964dcce0ee2SMatthew G. Knepley         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j];
965494e7359SMatthew G. Knepley       }
966494e7359SMatthew G. Knepley     }
967494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
968494e7359SMatthew G. Knepley     break;
969494e7359SMatthew G. Knepley   case 3:
970dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr);
971dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
972dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
973dcce0ee2SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
974dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
975dcce0ee2SMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
976dcce0ee2SMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
977dcce0ee2SMatthew 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);
978dcce0ee2SMatthew G. Knepley           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k];
979494e7359SMatthew G. Knepley         }
980494e7359SMatthew G. Knepley       }
981494e7359SMatthew G. Knepley     }
982494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
983494e7359SMatthew G. Knepley     break;
984494e7359SMatthew G. Knepley   default:
985494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
986494e7359SMatthew G. Knepley   }
98721454ff5SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
9882f5fb066SToby Isaac   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
989dcce0ee2SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
990d9bac1caSLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr);
991494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
992494e7359SMatthew G. Knepley }
993494e7359SMatthew G. Knepley 
994f5f57ec0SBarry Smith /*@
995b3c0f97bSTom Klotz   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
996b3c0f97bSTom Klotz 
997b3c0f97bSTom Klotz   Not Collective
998b3c0f97bSTom Klotz 
999b3c0f97bSTom Klotz   Input Arguments:
1000b3c0f97bSTom Klotz + dim   - The cell dimension
1001b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l
1002b3c0f97bSTom Klotz . a     - left end of interval (often-1)
1003b3c0f97bSTom Klotz - b     - right end of interval (often +1)
1004b3c0f97bSTom Klotz 
1005b3c0f97bSTom Klotz   Output Argument:
1006b3c0f97bSTom Klotz . q - A PetscQuadrature object
1007b3c0f97bSTom Klotz 
1008b3c0f97bSTom Klotz   Level: intermediate
1009b3c0f97bSTom Klotz 
1010b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature()
1011b3c0f97bSTom Klotz @*/
1012b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
1013b3c0f97bSTom Klotz {
1014b3c0f97bSTom Klotz   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
1015b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;                  /* Half-width of the integration interval */
1016b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;                  /* Center of the integration interval */
1017b3c0f97bSTom Klotz   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
1018d84b4d08SMatthew G. Knepley   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
1019b3c0f97bSTom Klotz   PetscReal       wk    = 0.5*PETSC_PI;              /* Quadrature weight at x_k */
1020b3c0f97bSTom Klotz   PetscReal      *x, *w;
1021b3c0f97bSTom Klotz   PetscInt        K, k, npoints;
1022b3c0f97bSTom Klotz   PetscErrorCode  ierr;
1023b3c0f97bSTom Klotz 
1024b3c0f97bSTom Klotz   PetscFunctionBegin;
1025b3c0f97bSTom Klotz   if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
1026b3c0f97bSTom Klotz   if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
1027b3c0f97bSTom Klotz   /* Find K such that the weights are < 32 digits of precision */
1028b3c0f97bSTom Klotz   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
10299add2064SThomas Klotz     wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
1030b3c0f97bSTom Klotz   }
1031b3c0f97bSTom Klotz   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
1032b3c0f97bSTom Klotz   ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr);
1033b3c0f97bSTom Klotz   npoints = 2*K-1;
1034b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
1035b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
1036b3c0f97bSTom Klotz   /* Center term */
1037b3c0f97bSTom Klotz   x[0] = beta;
1038b3c0f97bSTom Klotz   w[0] = 0.5*alpha*PETSC_PI;
1039b3c0f97bSTom Klotz   for (k = 1; k < K; ++k) {
10409add2064SThomas Klotz     wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
10411118d4bcSLisandro Dalcin     xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h));
1042b3c0f97bSTom Klotz     x[2*k-1] = -alpha*xk+beta;
1043b3c0f97bSTom Klotz     w[2*k-1] = wk;
1044b3c0f97bSTom Klotz     x[2*k+0] =  alpha*xk+beta;
1045b3c0f97bSTom Klotz     w[2*k+0] = wk;
1046b3c0f97bSTom Klotz   }
1047a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr);
1048b3c0f97bSTom Klotz   PetscFunctionReturn(0);
1049b3c0f97bSTom Klotz }
1050b3c0f97bSTom Klotz 
1051b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1052b3c0f97bSTom Klotz {
1053b3c0f97bSTom Klotz   const PetscInt  p     = 16;        /* Digits of precision in the evaluation */
1054b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;  /* Half-width of the integration interval */
1055b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;  /* Center of the integration interval */
1056b3c0f97bSTom Klotz   PetscReal       h     = 1.0;       /* Step size, length between x_k */
1057b3c0f97bSTom Klotz   PetscInt        l     = 0;         /* Level of refinement, h = 2^{-l} */
1058b3c0f97bSTom Klotz   PetscReal       osum  = 0.0;       /* Integral on last level */
1059b3c0f97bSTom Klotz   PetscReal       psum  = 0.0;       /* Integral on the level before the last level */
1060b3c0f97bSTom Klotz   PetscReal       sum;               /* Integral on current level */
1061446c295cSMatthew G. Knepley   PetscReal       yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1062b3c0f97bSTom Klotz   PetscReal       lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1063b3c0f97bSTom Klotz   PetscReal       wk;                /* Quadrature weight at x_k */
1064b3c0f97bSTom Klotz   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
1065b3c0f97bSTom Klotz   PetscInt        d;                 /* Digits of precision in the integral */
1066b3c0f97bSTom Klotz 
1067b3c0f97bSTom Klotz   PetscFunctionBegin;
1068b3c0f97bSTom Klotz   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1069b3c0f97bSTom Klotz   /* Center term */
1070b3c0f97bSTom Klotz   func(beta, &lval);
1071b3c0f97bSTom Klotz   sum = 0.5*alpha*PETSC_PI*lval;
1072b3c0f97bSTom Klotz   /* */
1073b3c0f97bSTom Klotz   do {
1074b3c0f97bSTom Klotz     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
1075b3c0f97bSTom Klotz     PetscInt  k = 1;
1076b3c0f97bSTom Klotz 
1077b3c0f97bSTom Klotz     ++l;
1078b3c0f97bSTom Klotz     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1079b3c0f97bSTom Klotz     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1080b3c0f97bSTom Klotz     psum = osum;
1081b3c0f97bSTom Klotz     osum = sum;
1082b3c0f97bSTom Klotz     h   *= 0.5;
1083b3c0f97bSTom Klotz     sum *= 0.5;
1084b3c0f97bSTom Klotz     do {
10859add2064SThomas Klotz       wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1086446c295cSMatthew G. Knepley       yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1087446c295cSMatthew G. Knepley       lx = -alpha*(1.0 - yk)+beta;
1088446c295cSMatthew G. Knepley       rx =  alpha*(1.0 - yk)+beta;
1089b3c0f97bSTom Klotz       func(lx, &lval);
1090b3c0f97bSTom Klotz       func(rx, &rval);
1091b3c0f97bSTom Klotz       lterm   = alpha*wk*lval;
1092b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
1093b3c0f97bSTom Klotz       sum    += lterm;
1094b3c0f97bSTom Klotz       rterm   = alpha*wk*rval;
1095b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
1096b3c0f97bSTom Klotz       sum    += rterm;
1097b3c0f97bSTom Klotz       ++k;
1098b3c0f97bSTom Klotz       /* Only need to evaluate every other point on refined levels */
1099b3c0f97bSTom Klotz       if (l != 1) ++k;
11009add2064SThomas Klotz     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
1101b3c0f97bSTom Klotz 
1102b3c0f97bSTom Klotz     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
1103b3c0f97bSTom Klotz     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
1104b3c0f97bSTom Klotz     d3 = PetscLog10Real(maxTerm) - p;
110509d48545SBarry Smith     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
110609d48545SBarry Smith     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
1107b3c0f97bSTom Klotz     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
11089add2064SThomas Klotz   } while (d < digits && l < 12);
1109b3c0f97bSTom Klotz   *sol = sum;
1110e510cb1fSThomas Klotz 
1111b3c0f97bSTom Klotz   PetscFunctionReturn(0);
1112b3c0f97bSTom Klotz }
1113b3c0f97bSTom Klotz 
1114497880caSRichard Tran Mills #if defined(PETSC_HAVE_MPFR)
111529f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
111629f144ccSMatthew G. Knepley {
1117e510cb1fSThomas Klotz   const PetscInt  safetyFactor = 2;  /* Calculate abcissa until 2*p digits */
111829f144ccSMatthew G. Knepley   PetscInt        l            = 0;  /* Level of refinement, h = 2^{-l} */
111929f144ccSMatthew G. Knepley   mpfr_t          alpha;             /* Half-width of the integration interval */
112029f144ccSMatthew G. Knepley   mpfr_t          beta;              /* Center of the integration interval */
112129f144ccSMatthew G. Knepley   mpfr_t          h;                 /* Step size, length between x_k */
112229f144ccSMatthew G. Knepley   mpfr_t          osum;              /* Integral on last level */
112329f144ccSMatthew G. Knepley   mpfr_t          psum;              /* Integral on the level before the last level */
112429f144ccSMatthew G. Knepley   mpfr_t          sum;               /* Integral on current level */
112529f144ccSMatthew G. Knepley   mpfr_t          yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
112629f144ccSMatthew G. Knepley   mpfr_t          lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
112729f144ccSMatthew G. Knepley   mpfr_t          wk;                /* Quadrature weight at x_k */
112829f144ccSMatthew G. Knepley   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
112929f144ccSMatthew G. Knepley   PetscInt        d;                 /* Digits of precision in the integral */
113029f144ccSMatthew G. Knepley   mpfr_t          pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
113129f144ccSMatthew G. Knepley 
113229f144ccSMatthew G. Knepley   PetscFunctionBegin;
113329f144ccSMatthew G. Knepley   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
113429f144ccSMatthew G. Knepley   /* Create high precision storage */
1135c9f744b5SMatthew 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);
113629f144ccSMatthew G. Knepley   /* Initialization */
113729f144ccSMatthew G. Knepley   mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
113829f144ccSMatthew G. Knepley   mpfr_set_d(beta,  0.5*(b+a), MPFR_RNDN);
113929f144ccSMatthew G. Knepley   mpfr_set_d(osum,  0.0,       MPFR_RNDN);
114029f144ccSMatthew G. Knepley   mpfr_set_d(psum,  0.0,       MPFR_RNDN);
114129f144ccSMatthew G. Knepley   mpfr_set_d(h,     1.0,       MPFR_RNDN);
114229f144ccSMatthew G. Knepley   mpfr_const_pi(pi2, MPFR_RNDN);
114329f144ccSMatthew G. Knepley   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
114429f144ccSMatthew G. Knepley   /* Center term */
114529f144ccSMatthew G. Knepley   func(0.5*(b+a), &lval);
114629f144ccSMatthew G. Knepley   mpfr_set(sum, pi2, MPFR_RNDN);
114729f144ccSMatthew G. Knepley   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
114829f144ccSMatthew G. Knepley   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
114929f144ccSMatthew G. Knepley   /* */
115029f144ccSMatthew G. Knepley   do {
115129f144ccSMatthew G. Knepley     PetscReal d1, d2, d3, d4;
115229f144ccSMatthew G. Knepley     PetscInt  k = 1;
115329f144ccSMatthew G. Knepley 
115429f144ccSMatthew G. Knepley     ++l;
115529f144ccSMatthew G. Knepley     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
115629f144ccSMatthew G. Knepley     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
115729f144ccSMatthew G. Knepley     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
115829f144ccSMatthew G. Knepley     mpfr_set(psum, osum, MPFR_RNDN);
115929f144ccSMatthew G. Knepley     mpfr_set(osum,  sum, MPFR_RNDN);
116029f144ccSMatthew G. Knepley     mpfr_mul_d(h,   h,   0.5, MPFR_RNDN);
116129f144ccSMatthew G. Knepley     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
116229f144ccSMatthew G. Knepley     do {
116329f144ccSMatthew G. Knepley       mpfr_set_si(kh, k, MPFR_RNDN);
116429f144ccSMatthew G. Knepley       mpfr_mul(kh, kh, h, MPFR_RNDN);
116529f144ccSMatthew G. Knepley       /* Weight */
116629f144ccSMatthew G. Knepley       mpfr_set(wk, h, MPFR_RNDN);
116729f144ccSMatthew G. Knepley       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
116829f144ccSMatthew G. Knepley       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
116929f144ccSMatthew G. Knepley       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
117029f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
117129f144ccSMatthew G. Knepley       mpfr_sqr(tmp, tmp, MPFR_RNDN);
117229f144ccSMatthew G. Knepley       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
117329f144ccSMatthew G. Knepley       mpfr_div(wk, wk, tmp, MPFR_RNDN);
117429f144ccSMatthew G. Knepley       /* Abscissa */
117529f144ccSMatthew G. Knepley       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
117629f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
117729f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
117829f144ccSMatthew G. Knepley       mpfr_exp(tmp, msinh, MPFR_RNDN);
117929f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
118029f144ccSMatthew G. Knepley       /* Quadrature points */
118129f144ccSMatthew G. Knepley       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
118229f144ccSMatthew G. Knepley       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
118329f144ccSMatthew G. Knepley       mpfr_add(lx, lx, beta, MPFR_RNDU);
118429f144ccSMatthew G. Knepley       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
118529f144ccSMatthew G. Knepley       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
118629f144ccSMatthew G. Knepley       mpfr_add(rx, rx, beta, MPFR_RNDD);
118729f144ccSMatthew G. Knepley       /* Evaluation */
118829f144ccSMatthew G. Knepley       func(mpfr_get_d(lx, MPFR_RNDU), &lval);
118929f144ccSMatthew G. Knepley       func(mpfr_get_d(rx, MPFR_RNDD), &rval);
119029f144ccSMatthew G. Knepley       /* Update */
119129f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
119229f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
119329f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
119429f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
119529f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
119629f144ccSMatthew G. Knepley       mpfr_set(curTerm, tmp, MPFR_RNDN);
119729f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
119829f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
119929f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
120029f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
120129f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
120229f144ccSMatthew G. Knepley       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
120329f144ccSMatthew G. Knepley       ++k;
120429f144ccSMatthew G. Knepley       /* Only need to evaluate every other point on refined levels */
120529f144ccSMatthew G. Knepley       if (l != 1) ++k;
120629f144ccSMatthew G. Knepley       mpfr_log10(tmp, wk, MPFR_RNDN);
120729f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
1208c9f744b5SMatthew G. Knepley     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
120929f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
121029f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
121129f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
121229f144ccSMatthew G. Knepley     d1 = mpfr_get_d(tmp, MPFR_RNDN);
121329f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
121429f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
121529f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
121629f144ccSMatthew G. Knepley     d2 = mpfr_get_d(tmp, MPFR_RNDN);
121729f144ccSMatthew G. Knepley     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
1218c9f744b5SMatthew G. Knepley     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
121929f144ccSMatthew G. Knepley     mpfr_log10(tmp, curTerm, MPFR_RNDN);
122029f144ccSMatthew G. Knepley     d4 = mpfr_get_d(tmp, MPFR_RNDN);
122129f144ccSMatthew G. Knepley     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1222b0649871SThomas Klotz   } while (d < digits && l < 8);
122329f144ccSMatthew G. Knepley   *sol = mpfr_get_d(sum, MPFR_RNDN);
122429f144ccSMatthew G. Knepley   /* Cleanup */
122529f144ccSMatthew G. Knepley   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
122629f144ccSMatthew G. Knepley   PetscFunctionReturn(0);
122729f144ccSMatthew G. Knepley }
1228d525116cSMatthew G. Knepley #else
1229fbfcfee5SBarry Smith 
1230d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1231d525116cSMatthew G. Knepley {
1232d525116cSMatthew G. Knepley   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
1233d525116cSMatthew G. Knepley }
123429f144ccSMatthew G. Knepley #endif
123529f144ccSMatthew G. Knepley 
1236194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
1237194825f6SJed Brown  * A in column-major format
1238194825f6SJed Brown  * Ainv in row-major format
1239194825f6SJed Brown  * tau has length m
1240194825f6SJed Brown  * worksize must be >= max(1,n)
1241194825f6SJed Brown  */
1242194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1243194825f6SJed Brown {
1244194825f6SJed Brown   PetscErrorCode ierr;
1245194825f6SJed Brown   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
1246194825f6SJed Brown   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
1247194825f6SJed Brown 
1248194825f6SJed Brown   PetscFunctionBegin;
1249194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
1250194825f6SJed Brown   {
1251194825f6SJed Brown     PetscInt i,j;
1252dcca6d9dSJed Brown     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
1253194825f6SJed Brown     for (j=0; j<n; j++) {
1254194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1255194825f6SJed Brown     }
1256194825f6SJed Brown     mstride = m;
1257194825f6SJed Brown   }
1258194825f6SJed Brown #else
1259194825f6SJed Brown   A = A_in;
1260194825f6SJed Brown   Ainv = Ainv_out;
1261194825f6SJed Brown #endif
1262194825f6SJed Brown 
1263194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
1264194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
1265194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
1266194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
1267194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1268001a771dSBarry Smith   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1269194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1270194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1271194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1272194825f6SJed Brown 
1273194825f6SJed Brown   /* Extract an explicit representation of Q */
1274194825f6SJed Brown   Q = Ainv;
1275194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
1276194825f6SJed Brown   K = N;                        /* full rank */
1277c964aadfSJose E. Roman   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1278194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1279194825f6SJed Brown 
1280194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1281194825f6SJed Brown   Alpha = 1.0;
1282194825f6SJed Brown   ldb = lda;
1283001a771dSBarry Smith   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1284194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
1285194825f6SJed Brown 
1286194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
1287194825f6SJed Brown   {
1288194825f6SJed Brown     PetscInt i;
1289194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1290194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
1291194825f6SJed Brown   }
1292194825f6SJed Brown #endif
1293194825f6SJed Brown   PetscFunctionReturn(0);
1294194825f6SJed Brown }
1295194825f6SJed Brown 
1296194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1297194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1298194825f6SJed Brown {
1299194825f6SJed Brown   PetscErrorCode ierr;
1300194825f6SJed Brown   PetscReal      *Bv;
1301194825f6SJed Brown   PetscInt       i,j;
1302194825f6SJed Brown 
1303194825f6SJed Brown   PetscFunctionBegin;
1304785e854fSJed Brown   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
1305194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
1306194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
1307194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1308194825f6SJed Brown   for (i=0; i<ninterval; i++) {
1309194825f6SJed Brown     for (j=0; j<ndegree; j++) {
1310194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1311194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1312194825f6SJed Brown     }
1313194825f6SJed Brown   }
1314194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
1315194825f6SJed Brown   PetscFunctionReturn(0);
1316194825f6SJed Brown }
1317194825f6SJed Brown 
1318194825f6SJed Brown /*@
1319194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
1320194825f6SJed Brown 
1321194825f6SJed Brown    Not Collective
1322194825f6SJed Brown 
1323194825f6SJed Brown    Input Arguments:
1324194825f6SJed Brown +  degree - degree of reconstruction polynomial
1325194825f6SJed Brown .  nsource - number of source intervals
1326194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1327194825f6SJed Brown .  ntarget - number of target intervals
1328194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
1329194825f6SJed Brown 
1330194825f6SJed Brown    Output Arguments:
1331194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
1332194825f6SJed Brown 
1333194825f6SJed Brown    Level: advanced
1334194825f6SJed Brown 
1335194825f6SJed Brown .seealso: PetscDTLegendreEval()
1336194825f6SJed Brown @*/
1337194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1338194825f6SJed Brown {
1339194825f6SJed Brown   PetscErrorCode ierr;
1340194825f6SJed Brown   PetscInt       i,j,k,*bdegrees,worksize;
1341194825f6SJed Brown   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1342194825f6SJed Brown   PetscScalar    *tau,*work;
1343194825f6SJed Brown 
1344194825f6SJed Brown   PetscFunctionBegin;
1345194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
1346194825f6SJed Brown   PetscValidRealPointer(targetx,5);
1347194825f6SJed Brown   PetscValidRealPointer(R,6);
1348194825f6SJed 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);
1349194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
1350194825f6SJed Brown   for (i=0; i<nsource; i++) {
135157622a8eSBarry 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]);
1352194825f6SJed Brown   }
1353194825f6SJed Brown   for (i=0; i<ntarget; i++) {
135457622a8eSBarry 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]);
1355194825f6SJed Brown   }
1356194825f6SJed Brown #endif
1357194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
1358194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1359194825f6SJed Brown   center = (xmin + xmax)/2;
1360194825f6SJed Brown   hscale = (xmax - xmin)/2;
1361194825f6SJed Brown   worksize = nsource;
1362dcca6d9dSJed Brown   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
1363dcca6d9dSJed Brown   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
1364194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1365194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1366194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
1367194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
1368194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1369194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
1370194825f6SJed Brown   for (i=0; i<ntarget; i++) {
1371194825f6SJed Brown     PetscReal rowsum = 0;
1372194825f6SJed Brown     for (j=0; j<nsource; j++) {
1373194825f6SJed Brown       PetscReal sum = 0;
1374194825f6SJed Brown       for (k=0; k<degree+1; k++) {
1375194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1376194825f6SJed Brown       }
1377194825f6SJed Brown       R[i*nsource+j] = sum;
1378194825f6SJed Brown       rowsum += sum;
1379194825f6SJed Brown     }
1380194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1381194825f6SJed Brown   }
1382194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
1383194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
1384194825f6SJed Brown   PetscFunctionReturn(0);
1385194825f6SJed Brown }
1386*916e780bShannah_mairs 
1387*916e780bShannah_mairs /*@C
1388*916e780bShannah_mairs    PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points
1389*916e780bShannah_mairs 
1390*916e780bShannah_mairs    Not Collective
1391*916e780bShannah_mairs 
1392*916e780bShannah_mairs    Input Parameter:
1393*916e780bShannah_mairs +  n - the number of GLL nodes
1394*916e780bShannah_mairs .  nodes - the GLL nodes
1395*916e780bShannah_mairs .  weights - the GLL weights
1396*916e780bShannah_mairs .  f - the function values at the nodes
1397*916e780bShannah_mairs 
1398*916e780bShannah_mairs    Output Parameter:
1399*916e780bShannah_mairs .  in - the value of the integral
1400*916e780bShannah_mairs 
1401*916e780bShannah_mairs    Level: beginner
1402*916e780bShannah_mairs 
1403*916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature()
1404*916e780bShannah_mairs 
1405*916e780bShannah_mairs @*/
1406*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in)
1407*916e780bShannah_mairs {
1408*916e780bShannah_mairs   PetscInt          i;
1409*916e780bShannah_mairs 
1410*916e780bShannah_mairs   PetscFunctionBegin;
1411*916e780bShannah_mairs   *in = 0.;
1412*916e780bShannah_mairs   for (i=0; i<n; i++) {
1413*916e780bShannah_mairs     *in += f[i]*f[i]*weights[i];
1414*916e780bShannah_mairs   }
1415*916e780bShannah_mairs   PetscFunctionReturn(0);
1416*916e780bShannah_mairs }
1417*916e780bShannah_mairs 
1418*916e780bShannah_mairs static void gllqAndLEvaluation(PetscInt n,PetscReal x,PetscReal *q,PetscReal *qp,PetscReal *Ln)
1419*916e780bShannah_mairs /*
1420*916e780bShannah_mairs   Compute the polynomial q(x) = L_{N+1}(x) - L_{n-1}(x) and its derivative in
1421*916e780bShannah_mairs   addition to L_N(x) as these are needed for computing the GLL points via Newton's method.
1422*916e780bShannah_mairs   Reference: "Implementing Spectral Methods for Partial Differential Equations: Algorithms
1423*916e780bShannah_mairs   for Scientists and Engineers" by David A. Kopriva.
1424*916e780bShannah_mairs */
1425*916e780bShannah_mairs {
1426*916e780bShannah_mairs   PetscInt k;
1427*916e780bShannah_mairs 
1428*916e780bShannah_mairs   PetscReal Lnp;
1429*916e780bShannah_mairs   PetscReal Lnp1, Lnp1p;
1430*916e780bShannah_mairs   PetscReal Lnm1, Lnm1p;
1431*916e780bShannah_mairs   PetscReal Lnm2, Lnm2p;
1432*916e780bShannah_mairs 
1433*916e780bShannah_mairs   Lnm1  = 1.0;
1434*916e780bShannah_mairs   *Ln   = x;
1435*916e780bShannah_mairs   Lnm1p = 0.0;
1436*916e780bShannah_mairs   Lnp   = 1.0;
1437*916e780bShannah_mairs 
1438*916e780bShannah_mairs   for (k=2; k<=n; ++k) {
1439*916e780bShannah_mairs     Lnm2  = Lnm1;
1440*916e780bShannah_mairs     Lnm1  = *Ln;
1441*916e780bShannah_mairs     Lnm2p = Lnm1p;
1442*916e780bShannah_mairs     Lnm1p = Lnp;
1443*916e780bShannah_mairs     *Ln   = (2.*((PetscReal)k)-1.)/(1.0*((PetscReal)k))*x*Lnm1 - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm2;
1444*916e780bShannah_mairs     Lnp   = Lnm2p + (2.0*((PetscReal)k)-1.)*Lnm1;
1445*916e780bShannah_mairs   }
1446*916e780bShannah_mairs   k     = n+1;
1447*916e780bShannah_mairs   Lnp1  = (2.*((PetscReal)k)-1.)/(((PetscReal)k))*x*(*Ln) - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm1;
1448*916e780bShannah_mairs   Lnp1p = Lnm1p + (2.0*((PetscReal)k)-1.)*(*Ln);
1449*916e780bShannah_mairs   *q    = Lnp1 - Lnm1;
1450*916e780bShannah_mairs   *qp   = Lnp1p - Lnm1p;
1451*916e780bShannah_mairs }
1452*916e780bShannah_mairs 
1453*916e780bShannah_mairs /*@C
1454*916e780bShannah_mairs    PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element
1455*916e780bShannah_mairs 
1456*916e780bShannah_mairs    Not Collective
1457*916e780bShannah_mairs 
1458*916e780bShannah_mairs    Input Parameter:
1459*916e780bShannah_mairs +  n - the number of GLL nodes
1460*916e780bShannah_mairs .  nodes - the GLL nodes
1461*916e780bShannah_mairs .  weights - the GLL weights
1462*916e780bShannah_mairs 
1463*916e780bShannah_mairs    Output Parameter:
1464*916e780bShannah_mairs .  A - the stiffness element
1465*916e780bShannah_mairs 
1466*916e780bShannah_mairs    Level: beginner
1467*916e780bShannah_mairs 
1468*916e780bShannah_mairs    Notes:
1469*916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy()
1470*916e780bShannah_mairs 
1471*916e780bShannah_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)
1472*916e780bShannah_mairs 
1473*916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1474*916e780bShannah_mairs 
1475*916e780bShannah_mairs @*/
1476*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1477*916e780bShannah_mairs {
1478*916e780bShannah_mairs   PetscReal        **A;
1479*916e780bShannah_mairs   PetscErrorCode  ierr;
1480*916e780bShannah_mairs   const PetscReal  *gllnodes = nodes;
1481*916e780bShannah_mairs   const PetscInt   p = n-1;
1482*916e780bShannah_mairs   PetscReal        z0,z1,z2 = -1,x,Lpj,Lpr;
1483*916e780bShannah_mairs   PetscInt         i,j,nn,r;
1484*916e780bShannah_mairs 
1485*916e780bShannah_mairs   PetscFunctionBegin;
1486*916e780bShannah_mairs   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1487*916e780bShannah_mairs   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1488*916e780bShannah_mairs   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1489*916e780bShannah_mairs 
1490*916e780bShannah_mairs   for (j=1; j<p; j++) {
1491*916e780bShannah_mairs     x  = gllnodes[j];
1492*916e780bShannah_mairs     z0 = 1.;
1493*916e780bShannah_mairs     z1 = x;
1494*916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1495*916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1496*916e780bShannah_mairs       z0 = z1;
1497*916e780bShannah_mairs       z1 = z2;
1498*916e780bShannah_mairs     }
1499*916e780bShannah_mairs     Lpj=z2;
1500*916e780bShannah_mairs     for (r=1; r<p; r++) {
1501*916e780bShannah_mairs       if (r == j) {
1502*916e780bShannah_mairs         A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj);
1503*916e780bShannah_mairs       } else {
1504*916e780bShannah_mairs         x  = gllnodes[r];
1505*916e780bShannah_mairs         z0 = 1.;
1506*916e780bShannah_mairs         z1 = x;
1507*916e780bShannah_mairs         for (nn=1; nn<p; nn++) {
1508*916e780bShannah_mairs           z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1509*916e780bShannah_mairs           z0 = z1;
1510*916e780bShannah_mairs           z1 = z2;
1511*916e780bShannah_mairs         }
1512*916e780bShannah_mairs         Lpr     = z2;
1513*916e780bShannah_mairs         A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r]));
1514*916e780bShannah_mairs       }
1515*916e780bShannah_mairs     }
1516*916e780bShannah_mairs   }
1517*916e780bShannah_mairs   for (j=1; j<p+1; j++) {
1518*916e780bShannah_mairs     x  = gllnodes[j];
1519*916e780bShannah_mairs     z0 = 1.;
1520*916e780bShannah_mairs     z1 = x;
1521*916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1522*916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1523*916e780bShannah_mairs       z0 = z1;
1524*916e780bShannah_mairs       z1 = z2;
1525*916e780bShannah_mairs     }
1526*916e780bShannah_mairs     Lpj     = z2;
1527*916e780bShannah_mairs     A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j]));
1528*916e780bShannah_mairs     A[0][j] = A[j][0];
1529*916e780bShannah_mairs   }
1530*916e780bShannah_mairs   for (j=0; j<p; j++) {
1531*916e780bShannah_mairs     x  = gllnodes[j];
1532*916e780bShannah_mairs     z0 = 1.;
1533*916e780bShannah_mairs     z1 = x;
1534*916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1535*916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1536*916e780bShannah_mairs       z0 = z1;
1537*916e780bShannah_mairs       z1 = z2;
1538*916e780bShannah_mairs     }
1539*916e780bShannah_mairs     Lpj=z2;
1540*916e780bShannah_mairs 
1541*916e780bShannah_mairs     A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j]));
1542*916e780bShannah_mairs     A[j][p] = A[p][j];
1543*916e780bShannah_mairs   }
1544*916e780bShannah_mairs   A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.;
1545*916e780bShannah_mairs   A[p][p]=A[0][0];
1546*916e780bShannah_mairs   *AA = A;
1547*916e780bShannah_mairs   PetscFunctionReturn(0);
1548*916e780bShannah_mairs }
1549*916e780bShannah_mairs 
1550*916e780bShannah_mairs /*@C
1551*916e780bShannah_mairs    PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element
1552*916e780bShannah_mairs 
1553*916e780bShannah_mairs    Not Collective
1554*916e780bShannah_mairs 
1555*916e780bShannah_mairs    Input Parameter:
1556*916e780bShannah_mairs +  n - the number of GLL nodes
1557*916e780bShannah_mairs .  nodes - the GLL nodes
1558*916e780bShannah_mairs .  weights - the GLL weightss
1559*916e780bShannah_mairs -  A - the stiffness element
1560*916e780bShannah_mairs 
1561*916e780bShannah_mairs    Level: beginner
1562*916e780bShannah_mairs 
1563*916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate()
1564*916e780bShannah_mairs 
1565*916e780bShannah_mairs @*/
1566*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1567*916e780bShannah_mairs {
1568*916e780bShannah_mairs   PetscErrorCode ierr;
1569*916e780bShannah_mairs 
1570*916e780bShannah_mairs   PetscFunctionBegin;
1571*916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1572*916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1573*916e780bShannah_mairs   *AA  = NULL;
1574*916e780bShannah_mairs   PetscFunctionReturn(0);
1575*916e780bShannah_mairs }
1576*916e780bShannah_mairs 
1577*916e780bShannah_mairs /*@C
1578*916e780bShannah_mairs    PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element
1579*916e780bShannah_mairs 
1580*916e780bShannah_mairs    Not Collective
1581*916e780bShannah_mairs 
1582*916e780bShannah_mairs    Input Parameter:
1583*916e780bShannah_mairs +  n - the number of GLL nodes
1584*916e780bShannah_mairs .  nodes - the GLL nodes
1585*916e780bShannah_mairs .  weights - the GLL weights
1586*916e780bShannah_mairs 
1587*916e780bShannah_mairs    Output Parameter:
1588*916e780bShannah_mairs .  AA - the stiffness element
1589*916e780bShannah_mairs -  AAT - the transpose of AA (pass in NULL if you do not need this array)
1590*916e780bShannah_mairs 
1591*916e780bShannah_mairs    Level: beginner
1592*916e780bShannah_mairs 
1593*916e780bShannah_mairs    Notes:
1594*916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementGradientDestroy()
1595*916e780bShannah_mairs 
1596*916e780bShannah_mairs    You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1597*916e780bShannah_mairs 
1598*916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1599*916e780bShannah_mairs 
1600*916e780bShannah_mairs @*/
1601*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
1602*916e780bShannah_mairs {
1603*916e780bShannah_mairs   PetscReal        **A, **AT = NULL;
1604*916e780bShannah_mairs   PetscErrorCode  ierr;
1605*916e780bShannah_mairs   const PetscReal  *gllnodes = nodes;
1606*916e780bShannah_mairs   const PetscInt   p = n-1;
1607*916e780bShannah_mairs   PetscReal        q,qp,Li, Lj,d0;
1608*916e780bShannah_mairs   PetscInt         i,j;
1609*916e780bShannah_mairs 
1610*916e780bShannah_mairs   PetscFunctionBegin;
1611*916e780bShannah_mairs   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1612*916e780bShannah_mairs   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1613*916e780bShannah_mairs   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1614*916e780bShannah_mairs 
1615*916e780bShannah_mairs   if (AAT) {
1616*916e780bShannah_mairs     ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr);
1617*916e780bShannah_mairs     ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr);
1618*916e780bShannah_mairs     for (i=1; i<n; i++) AT[i] = AT[i-1]+n;
1619*916e780bShannah_mairs   }
1620*916e780bShannah_mairs 
1621*916e780bShannah_mairs   if (n==1) {A[0][0] = 0.;}
1622*916e780bShannah_mairs   d0 = (PetscReal)p*((PetscReal)p+1.)/4.;
1623*916e780bShannah_mairs   for  (i=0; i<n; i++) {
1624*916e780bShannah_mairs     for  (j=0; j<n; j++) {
1625*916e780bShannah_mairs       A[i][j] = 0.;
1626*916e780bShannah_mairs       gllqAndLEvaluation(p,gllnodes[i],&q,&qp,&Li);
1627*916e780bShannah_mairs       gllqAndLEvaluation(p,gllnodes[j],&q,&qp,&Lj);
1628*916e780bShannah_mairs       if (i!=j)             A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j]));
1629*916e780bShannah_mairs       if ((j==i) && (i==0)) A[i][j] = -d0;
1630*916e780bShannah_mairs       if (j==i && i==p)     A[i][j] = d0;
1631*916e780bShannah_mairs       if (AT) AT[j][i] = A[i][j];
1632*916e780bShannah_mairs     }
1633*916e780bShannah_mairs   }
1634*916e780bShannah_mairs   if (AAT) *AAT = AT;
1635*916e780bShannah_mairs   *AA  = A;
1636*916e780bShannah_mairs   PetscFunctionReturn(0);
1637*916e780bShannah_mairs }
1638*916e780bShannah_mairs 
1639*916e780bShannah_mairs /*@C
1640*916e780bShannah_mairs    PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate()
1641*916e780bShannah_mairs 
1642*916e780bShannah_mairs    Not Collective
1643*916e780bShannah_mairs 
1644*916e780bShannah_mairs    Input Parameter:
1645*916e780bShannah_mairs +  n - the number of GLL nodes
1646*916e780bShannah_mairs .  nodes - the GLL nodes
1647*916e780bShannah_mairs .  weights - the GLL weights
1648*916e780bShannah_mairs .  AA - the stiffness element
1649*916e780bShannah_mairs -  AAT - the transpose of the element
1650*916e780bShannah_mairs 
1651*916e780bShannah_mairs    Level: beginner
1652*916e780bShannah_mairs 
1653*916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate()
1654*916e780bShannah_mairs 
1655*916e780bShannah_mairs @*/
1656*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
1657*916e780bShannah_mairs {
1658*916e780bShannah_mairs   PetscErrorCode ierr;
1659*916e780bShannah_mairs 
1660*916e780bShannah_mairs   PetscFunctionBegin;
1661*916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1662*916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1663*916e780bShannah_mairs   *AA  = NULL;
1664*916e780bShannah_mairs   if (*AAT) {
1665*916e780bShannah_mairs     ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr);
1666*916e780bShannah_mairs     ierr = PetscFree(*AAT);CHKERRQ(ierr);
1667*916e780bShannah_mairs     *AAT  = NULL;
1668*916e780bShannah_mairs   }
1669*916e780bShannah_mairs   PetscFunctionReturn(0);
1670*916e780bShannah_mairs }
1671*916e780bShannah_mairs 
1672*916e780bShannah_mairs /*@C
1673*916e780bShannah_mairs    PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element
1674*916e780bShannah_mairs 
1675*916e780bShannah_mairs    Not Collective
1676*916e780bShannah_mairs 
1677*916e780bShannah_mairs    Input Parameter:
1678*916e780bShannah_mairs +  n - the number of GLL nodes
1679*916e780bShannah_mairs .  nodes - the GLL nodes
1680*916e780bShannah_mairs .  weights - the GLL weightss
1681*916e780bShannah_mairs 
1682*916e780bShannah_mairs    Output Parameter:
1683*916e780bShannah_mairs .  AA - the stiffness element
1684*916e780bShannah_mairs 
1685*916e780bShannah_mairs    Level: beginner
1686*916e780bShannah_mairs 
1687*916e780bShannah_mairs    Notes:
1688*916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy()
1689*916e780bShannah_mairs 
1690*916e780bShannah_mairs    This is the same as the Gradient operator multiplied by the diagonal mass matrix
1691*916e780bShannah_mairs 
1692*916e780bShannah_mairs    You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1693*916e780bShannah_mairs 
1694*916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy()
1695*916e780bShannah_mairs 
1696*916e780bShannah_mairs @*/
1697*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1698*916e780bShannah_mairs {
1699*916e780bShannah_mairs   PetscReal       **D;
1700*916e780bShannah_mairs   PetscErrorCode  ierr;
1701*916e780bShannah_mairs   const PetscReal  *gllweights = weights;
1702*916e780bShannah_mairs   const PetscInt   glln = n;
1703*916e780bShannah_mairs   PetscInt         i,j;
1704*916e780bShannah_mairs 
1705*916e780bShannah_mairs   PetscFunctionBegin;
1706*916e780bShannah_mairs   ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr);
1707*916e780bShannah_mairs   for (i=0; i<glln; i++){
1708*916e780bShannah_mairs     for (j=0; j<glln; j++) {
1709*916e780bShannah_mairs       D[i][j] = gllweights[i]*D[i][j];
1710*916e780bShannah_mairs     }
1711*916e780bShannah_mairs   }
1712*916e780bShannah_mairs   *AA = D;
1713*916e780bShannah_mairs   PetscFunctionReturn(0);
1714*916e780bShannah_mairs }
1715*916e780bShannah_mairs 
1716*916e780bShannah_mairs /*@C
1717*916e780bShannah_mairs    PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element
1718*916e780bShannah_mairs 
1719*916e780bShannah_mairs    Not Collective
1720*916e780bShannah_mairs 
1721*916e780bShannah_mairs    Input Parameter:
1722*916e780bShannah_mairs +  n - the number of GLL nodes
1723*916e780bShannah_mairs .  nodes - the GLL nodes
1724*916e780bShannah_mairs .  weights - the GLL weights
1725*916e780bShannah_mairs -  A - advection
1726*916e780bShannah_mairs 
1727*916e780bShannah_mairs    Level: beginner
1728*916e780bShannah_mairs 
1729*916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate()
1730*916e780bShannah_mairs 
1731*916e780bShannah_mairs @*/
1732*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1733*916e780bShannah_mairs {
1734*916e780bShannah_mairs   PetscErrorCode ierr;
1735*916e780bShannah_mairs 
1736*916e780bShannah_mairs   PetscFunctionBegin;
1737*916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1738*916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1739*916e780bShannah_mairs   *AA  = NULL;
1740*916e780bShannah_mairs   PetscFunctionReturn(0);
1741*916e780bShannah_mairs }
1742*916e780bShannah_mairs 
1743*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1744*916e780bShannah_mairs {
1745*916e780bShannah_mairs   PetscReal        **A;
1746*916e780bShannah_mairs   PetscErrorCode  ierr;
1747*916e780bShannah_mairs   const PetscReal  *gllweights = weights;
1748*916e780bShannah_mairs   const PetscInt   glln = n;
1749*916e780bShannah_mairs   PetscInt         i,j;
1750*916e780bShannah_mairs 
1751*916e780bShannah_mairs   PetscFunctionBegin;
1752*916e780bShannah_mairs   ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr);
1753*916e780bShannah_mairs   ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr);
1754*916e780bShannah_mairs   for (i=1; i<glln; i++) A[i] = A[i-1]+glln;
1755*916e780bShannah_mairs   if (glln==1) {A[0][0] = 0.;}
1756*916e780bShannah_mairs   for  (i=0; i<glln; i++) {
1757*916e780bShannah_mairs     for  (j=0; j<glln; j++) {
1758*916e780bShannah_mairs       A[i][j] = 0.;
1759*916e780bShannah_mairs       if (j==i)     A[i][j] = gllweights[i];
1760*916e780bShannah_mairs     }
1761*916e780bShannah_mairs   }
1762*916e780bShannah_mairs   *AA  = A;
1763*916e780bShannah_mairs   PetscFunctionReturn(0);
1764*916e780bShannah_mairs }
1765*916e780bShannah_mairs 
1766*916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1767*916e780bShannah_mairs {
1768*916e780bShannah_mairs   PetscErrorCode ierr;
1769*916e780bShannah_mairs 
1770*916e780bShannah_mairs   PetscFunctionBegin;
1771*916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1772*916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1773*916e780bShannah_mairs   *AA  = NULL;
1774*916e780bShannah_mairs   PetscFunctionReturn(0);
1775*916e780bShannah_mairs }
1776*916e780bShannah_mairs 
1777