xref: /petsc/src/dm/dt/interface/dt.c (revision e6a796c3f68c2706b50e6802ecb449eebddacf45)
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 
15*e6a796c3SToby Isaac static PetscBool GolubWelschCite       = PETSC_FALSE;
16*e6a796c3SToby Isaac const char       GolubWelschCitation[] = "@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 
25*e6a796c3SToby Isaac static PetscBool gaussQuadratureNewton = PETSC_FALSE;
26*e6a796c3SToby Isaac 
272cd22861SMatthew G. Knepley 
282cd22861SMatthew G. Knepley PetscClassId PETSCQUADRATURE_CLASSID = 0;
292cd22861SMatthew G. Knepley 
3040d8ff71SMatthew G. Knepley /*@
3140d8ff71SMatthew G. Knepley   PetscQuadratureCreate - Create a PetscQuadrature object
3240d8ff71SMatthew G. Knepley 
33d083f849SBarry Smith   Collective
3440d8ff71SMatthew G. Knepley 
3540d8ff71SMatthew G. Knepley   Input Parameter:
3640d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object
3740d8ff71SMatthew G. Knepley 
3840d8ff71SMatthew G. Knepley   Output Parameter:
3940d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
4040d8ff71SMatthew G. Knepley 
4140d8ff71SMatthew G. Knepley   Level: beginner
4240d8ff71SMatthew G. Knepley 
4340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
4440d8ff71SMatthew G. Knepley @*/
4521454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
4621454ff5SMatthew G. Knepley {
4721454ff5SMatthew G. Knepley   PetscErrorCode ierr;
4821454ff5SMatthew G. Knepley 
4921454ff5SMatthew G. Knepley   PetscFunctionBegin;
5021454ff5SMatthew G. Knepley   PetscValidPointer(q, 2);
512cd22861SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
522cd22861SMatthew G. Knepley   ierr = PetscHeaderCreate(*q,PETSCQUADRATURE_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
5321454ff5SMatthew G. Knepley   (*q)->dim       = -1;
54a6b92713SMatthew G. Knepley   (*q)->Nc        =  1;
55bcede257SMatthew G. Knepley   (*q)->order     = -1;
5621454ff5SMatthew G. Knepley   (*q)->numPoints = 0;
5721454ff5SMatthew G. Knepley   (*q)->points    = NULL;
5821454ff5SMatthew G. Knepley   (*q)->weights   = NULL;
5921454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
6021454ff5SMatthew G. Knepley }
6121454ff5SMatthew G. Knepley 
62c9638911SMatthew G. Knepley /*@
63c9638911SMatthew G. Knepley   PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object
64c9638911SMatthew G. Knepley 
65d083f849SBarry Smith   Collective on q
66c9638911SMatthew G. Knepley 
67c9638911SMatthew G. Knepley   Input Parameter:
68c9638911SMatthew G. Knepley . q  - The PetscQuadrature object
69c9638911SMatthew G. Knepley 
70c9638911SMatthew G. Knepley   Output Parameter:
71c9638911SMatthew G. Knepley . r  - The new PetscQuadrature object
72c9638911SMatthew G. Knepley 
73c9638911SMatthew G. Knepley   Level: beginner
74c9638911SMatthew G. Knepley 
75c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData()
76c9638911SMatthew G. Knepley @*/
77c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
78c9638911SMatthew G. Knepley {
79a6b92713SMatthew G. Knepley   PetscInt         order, dim, Nc, Nq;
80c9638911SMatthew G. Knepley   const PetscReal *points, *weights;
81c9638911SMatthew G. Knepley   PetscReal       *p, *w;
82c9638911SMatthew G. Knepley   PetscErrorCode   ierr;
83c9638911SMatthew G. Knepley 
84c9638911SMatthew G. Knepley   PetscFunctionBegin;
85c9638911SMatthew G. Knepley   PetscValidPointer(q, 2);
86c9638911SMatthew G. Knepley   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr);
87c9638911SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
88c9638911SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr);
89a6b92713SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
90c9638911SMatthew G. Knepley   ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr);
91f0a0bfafSMatthew G. Knepley   ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr);
92580bdb30SBarry Smith   ierr = PetscArraycpy(p, points, Nq*dim);CHKERRQ(ierr);
93580bdb30SBarry Smith   ierr = PetscArraycpy(w, weights, Nc * Nq);CHKERRQ(ierr);
94a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr);
95c9638911SMatthew G. Knepley   PetscFunctionReturn(0);
96c9638911SMatthew G. Knepley }
97c9638911SMatthew G. Knepley 
9840d8ff71SMatthew G. Knepley /*@
9940d8ff71SMatthew G. Knepley   PetscQuadratureDestroy - Destroys a PetscQuadrature object
10040d8ff71SMatthew G. Knepley 
101d083f849SBarry Smith   Collective on q
10240d8ff71SMatthew G. Knepley 
10340d8ff71SMatthew G. Knepley   Input Parameter:
10440d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
10540d8ff71SMatthew G. Knepley 
10640d8ff71SMatthew G. Knepley   Level: beginner
10740d8ff71SMatthew G. Knepley 
10840d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
10940d8ff71SMatthew G. Knepley @*/
110bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
111bfa639d9SMatthew G. Knepley {
112bfa639d9SMatthew G. Knepley   PetscErrorCode ierr;
113bfa639d9SMatthew G. Knepley 
114bfa639d9SMatthew G. Knepley   PetscFunctionBegin;
11521454ff5SMatthew G. Knepley   if (!*q) PetscFunctionReturn(0);
1162cd22861SMatthew G. Knepley   PetscValidHeaderSpecific((*q),PETSCQUADRATURE_CLASSID,1);
11721454ff5SMatthew G. Knepley   if (--((PetscObject)(*q))->refct > 0) {
11821454ff5SMatthew G. Knepley     *q = NULL;
11921454ff5SMatthew G. Knepley     PetscFunctionReturn(0);
12021454ff5SMatthew G. Knepley   }
12121454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
12221454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
12321454ff5SMatthew G. Knepley   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
12421454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
12521454ff5SMatthew G. Knepley }
12621454ff5SMatthew G. Knepley 
127bcede257SMatthew G. Knepley /*@
128a6b92713SMatthew G. Knepley   PetscQuadratureGetOrder - Return the order of the method
129bcede257SMatthew G. Knepley 
130bcede257SMatthew G. Knepley   Not collective
131bcede257SMatthew G. Knepley 
132bcede257SMatthew G. Knepley   Input Parameter:
133bcede257SMatthew G. Knepley . q - The PetscQuadrature object
134bcede257SMatthew G. Knepley 
135bcede257SMatthew G. Knepley   Output Parameter:
136bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
137bcede257SMatthew G. Knepley 
138bcede257SMatthew G. Knepley   Level: intermediate
139bcede257SMatthew G. Knepley 
140bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
141bcede257SMatthew G. Knepley @*/
142bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
143bcede257SMatthew G. Knepley {
144bcede257SMatthew G. Knepley   PetscFunctionBegin;
1452cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
146bcede257SMatthew G. Knepley   PetscValidPointer(order, 2);
147bcede257SMatthew G. Knepley   *order = q->order;
148bcede257SMatthew G. Knepley   PetscFunctionReturn(0);
149bcede257SMatthew G. Knepley }
150bcede257SMatthew G. Knepley 
151bcede257SMatthew G. Knepley /*@
152a6b92713SMatthew G. Knepley   PetscQuadratureSetOrder - Return the order of the method
153bcede257SMatthew G. Knepley 
154bcede257SMatthew G. Knepley   Not collective
155bcede257SMatthew G. Knepley 
156bcede257SMatthew G. Knepley   Input Parameters:
157bcede257SMatthew G. Knepley + q - The PetscQuadrature object
158bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
159bcede257SMatthew G. Knepley 
160bcede257SMatthew G. Knepley   Level: intermediate
161bcede257SMatthew G. Knepley 
162bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
163bcede257SMatthew G. Knepley @*/
164bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
165bcede257SMatthew G. Knepley {
166bcede257SMatthew G. Knepley   PetscFunctionBegin;
1672cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
168bcede257SMatthew G. Knepley   q->order = order;
169bcede257SMatthew G. Knepley   PetscFunctionReturn(0);
170bcede257SMatthew G. Knepley }
171bcede257SMatthew G. Knepley 
172a6b92713SMatthew G. Knepley /*@
173a6b92713SMatthew G. Knepley   PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated
174a6b92713SMatthew G. Knepley 
175a6b92713SMatthew G. Knepley   Not collective
176a6b92713SMatthew G. Knepley 
177a6b92713SMatthew G. Knepley   Input Parameter:
178a6b92713SMatthew G. Knepley . q - The PetscQuadrature object
179a6b92713SMatthew G. Knepley 
180a6b92713SMatthew G. Knepley   Output Parameter:
181a6b92713SMatthew G. Knepley . Nc - The number of components
182a6b92713SMatthew G. Knepley 
183a6b92713SMatthew G. Knepley   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
184a6b92713SMatthew G. Knepley 
185a6b92713SMatthew G. Knepley   Level: intermediate
186a6b92713SMatthew G. Knepley 
187a6b92713SMatthew G. Knepley .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
188a6b92713SMatthew G. Knepley @*/
189a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc)
190a6b92713SMatthew G. Knepley {
191a6b92713SMatthew G. Knepley   PetscFunctionBegin;
1922cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
193a6b92713SMatthew G. Knepley   PetscValidPointer(Nc, 2);
194a6b92713SMatthew G. Knepley   *Nc = q->Nc;
195a6b92713SMatthew G. Knepley   PetscFunctionReturn(0);
196a6b92713SMatthew G. Knepley }
197a6b92713SMatthew G. Knepley 
198a6b92713SMatthew G. Knepley /*@
199a6b92713SMatthew G. Knepley   PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated
200a6b92713SMatthew G. Knepley 
201a6b92713SMatthew G. Knepley   Not collective
202a6b92713SMatthew G. Knepley 
203a6b92713SMatthew G. Knepley   Input Parameters:
204a6b92713SMatthew G. Knepley + q  - The PetscQuadrature object
205a6b92713SMatthew G. Knepley - Nc - The number of components
206a6b92713SMatthew G. Knepley 
207a6b92713SMatthew G. Knepley   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
208a6b92713SMatthew G. Knepley 
209a6b92713SMatthew G. Knepley   Level: intermediate
210a6b92713SMatthew G. Knepley 
211a6b92713SMatthew G. Knepley .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
212a6b92713SMatthew G. Knepley @*/
213a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc)
214a6b92713SMatthew G. Knepley {
215a6b92713SMatthew G. Knepley   PetscFunctionBegin;
2162cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
217a6b92713SMatthew G. Knepley   q->Nc = Nc;
218a6b92713SMatthew G. Knepley   PetscFunctionReturn(0);
219a6b92713SMatthew G. Knepley }
220a6b92713SMatthew G. Knepley 
22140d8ff71SMatthew G. Knepley /*@C
22240d8ff71SMatthew G. Knepley   PetscQuadratureGetData - Returns the data defining the quadrature
22340d8ff71SMatthew G. Knepley 
22440d8ff71SMatthew G. Knepley   Not collective
22540d8ff71SMatthew G. Knepley 
22640d8ff71SMatthew G. Knepley   Input Parameter:
22740d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
22840d8ff71SMatthew G. Knepley 
22940d8ff71SMatthew G. Knepley   Output Parameters:
23040d8ff71SMatthew G. Knepley + dim - The spatial dimension
231805e7170SToby Isaac . Nc - The number of components
23240d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
23340d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
23440d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
23540d8ff71SMatthew G. Knepley 
23640d8ff71SMatthew G. Knepley   Level: intermediate
23740d8ff71SMatthew G. Knepley 
23895452b02SPatrick Sanan   Fortran Notes:
23995452b02SPatrick Sanan     From Fortran you must call PetscQuadratureRestoreData() when you are done with the data
2401fd49c25SBarry Smith 
24140d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
24240d8ff71SMatthew G. Knepley @*/
243a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
24421454ff5SMatthew G. Knepley {
24521454ff5SMatthew G. Knepley   PetscFunctionBegin;
2462cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
24721454ff5SMatthew G. Knepley   if (dim) {
24821454ff5SMatthew G. Knepley     PetscValidPointer(dim, 2);
24921454ff5SMatthew G. Knepley     *dim = q->dim;
25021454ff5SMatthew G. Knepley   }
251a6b92713SMatthew G. Knepley   if (Nc) {
252a6b92713SMatthew G. Knepley     PetscValidPointer(Nc, 3);
253a6b92713SMatthew G. Knepley     *Nc = q->Nc;
254a6b92713SMatthew G. Knepley   }
25521454ff5SMatthew G. Knepley   if (npoints) {
256a6b92713SMatthew G. Knepley     PetscValidPointer(npoints, 4);
25721454ff5SMatthew G. Knepley     *npoints = q->numPoints;
25821454ff5SMatthew G. Knepley   }
25921454ff5SMatthew G. Knepley   if (points) {
260a6b92713SMatthew G. Knepley     PetscValidPointer(points, 5);
26121454ff5SMatthew G. Knepley     *points = q->points;
26221454ff5SMatthew G. Knepley   }
26321454ff5SMatthew G. Knepley   if (weights) {
264a6b92713SMatthew G. Knepley     PetscValidPointer(weights, 6);
26521454ff5SMatthew G. Knepley     *weights = q->weights;
26621454ff5SMatthew G. Knepley   }
26721454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
26821454ff5SMatthew G. Knepley }
26921454ff5SMatthew G. Knepley 
270907761f8SToby Isaac static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[])
271907761f8SToby Isaac {
272907761f8SToby Isaac   PetscScalar    *Js, *Jinvs;
273907761f8SToby Isaac   PetscInt       i, j, k;
274907761f8SToby Isaac   PetscBLASInt   bm, bn, info;
275907761f8SToby Isaac   PetscErrorCode ierr;
276907761f8SToby Isaac 
277907761f8SToby Isaac   PetscFunctionBegin;
278907761f8SToby Isaac   ierr = PetscBLASIntCast(m, &bm);CHKERRQ(ierr);
279907761f8SToby Isaac   ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr);
280907761f8SToby Isaac #if defined(PETSC_USE_COMPLEX)
281907761f8SToby Isaac   ierr = PetscMalloc2(m*n, &Js, m*n, &Jinvs);CHKERRQ(ierr);
28228222859SToby Isaac   for (i = 0; i < m*n; i++) Js[i] = J[i];
283907761f8SToby Isaac #else
284907761f8SToby Isaac   Js = (PetscReal *) J;
285907761f8SToby Isaac   Jinvs = Jinv;
286907761f8SToby Isaac #endif
287907761f8SToby Isaac   if (m == n) {
288907761f8SToby Isaac     PetscBLASInt *pivots;
289907761f8SToby Isaac     PetscScalar *W;
290907761f8SToby Isaac 
291907761f8SToby Isaac     ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr);
292907761f8SToby Isaac 
293907761f8SToby Isaac     ierr = PetscArraycpy(Jinvs, Js, m * m);CHKERRQ(ierr);
294907761f8SToby Isaac     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info));
295907761f8SToby Isaac     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
296907761f8SToby Isaac     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info));
297907761f8SToby Isaac     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
298907761f8SToby Isaac     ierr = PetscFree2(pivots, W);CHKERRQ(ierr);
299907761f8SToby Isaac   } else if (m < n) {
300907761f8SToby Isaac     PetscScalar *JJT;
301907761f8SToby Isaac     PetscBLASInt *pivots;
302907761f8SToby Isaac     PetscScalar *W;
303907761f8SToby Isaac 
304907761f8SToby Isaac     ierr = PetscMalloc1(m*m, &JJT);CHKERRQ(ierr);
305907761f8SToby Isaac     ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr);
306907761f8SToby Isaac     for (i = 0; i < m; i++) {
307907761f8SToby Isaac       for (j = 0; j < m; j++) {
308907761f8SToby Isaac         PetscScalar val = 0.;
309907761f8SToby Isaac 
310907761f8SToby Isaac         for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k];
311907761f8SToby Isaac         JJT[i * m + j] = val;
312907761f8SToby Isaac       }
313907761f8SToby Isaac     }
314907761f8SToby Isaac 
315907761f8SToby Isaac     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info));
316907761f8SToby Isaac     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
317907761f8SToby Isaac     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info));
318907761f8SToby Isaac     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
319907761f8SToby Isaac     for (i = 0; i < n; i++) {
320907761f8SToby Isaac       for (j = 0; j < m; j++) {
321907761f8SToby Isaac         PetscScalar val = 0.;
322907761f8SToby Isaac 
323907761f8SToby Isaac         for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j];
324907761f8SToby Isaac         Jinvs[i * m + j] = val;
325907761f8SToby Isaac       }
326907761f8SToby Isaac     }
327907761f8SToby Isaac     ierr = PetscFree2(pivots, W);CHKERRQ(ierr);
328907761f8SToby Isaac     ierr = PetscFree(JJT);CHKERRQ(ierr);
329907761f8SToby Isaac   } else {
330907761f8SToby Isaac     PetscScalar *JTJ;
331907761f8SToby Isaac     PetscBLASInt *pivots;
332907761f8SToby Isaac     PetscScalar *W;
333907761f8SToby Isaac 
334907761f8SToby Isaac     ierr = PetscMalloc1(n*n, &JTJ);CHKERRQ(ierr);
335907761f8SToby Isaac     ierr = PetscMalloc2(n, &pivots, n, &W);CHKERRQ(ierr);
336907761f8SToby Isaac     for (i = 0; i < n; i++) {
337907761f8SToby Isaac       for (j = 0; j < n; j++) {
338907761f8SToby Isaac         PetscScalar val = 0.;
339907761f8SToby Isaac 
340907761f8SToby Isaac         for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j];
341907761f8SToby Isaac         JTJ[i * n + j] = val;
342907761f8SToby Isaac       }
343907761f8SToby Isaac     }
344907761f8SToby Isaac 
345907761f8SToby Isaac     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bm, pivots, &info));
346907761f8SToby Isaac     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
347907761f8SToby Isaac     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info));
348907761f8SToby Isaac     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
349907761f8SToby Isaac     for (i = 0; i < n; i++) {
350907761f8SToby Isaac       for (j = 0; j < m; j++) {
351907761f8SToby Isaac         PetscScalar val = 0.;
352907761f8SToby Isaac 
353907761f8SToby Isaac         for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k];
354907761f8SToby Isaac         Jinvs[i * m + j] = val;
355907761f8SToby Isaac       }
356907761f8SToby Isaac     }
357907761f8SToby Isaac     ierr = PetscFree2(pivots, W);CHKERRQ(ierr);
358907761f8SToby Isaac     ierr = PetscFree(JTJ);CHKERRQ(ierr);
359907761f8SToby Isaac   }
360907761f8SToby Isaac #if defined(PETSC_USE_COMPLEX)
36128222859SToby Isaac   for (i = 0; i < m*n; i++) Jinv[i] = PetscRealPart(Jinvs[i]);
362907761f8SToby Isaac   ierr = PetscFree2(Js, Jinvs);CHKERRQ(ierr);
363907761f8SToby Isaac #endif
364907761f8SToby Isaac   PetscFunctionReturn(0);
365907761f8SToby Isaac }
366907761f8SToby Isaac 
367907761f8SToby Isaac /*@
368907761f8SToby Isaac    PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation.
369907761f8SToby Isaac 
370907761f8SToby Isaac    Collecive on PetscQuadrature
371907761f8SToby Isaac 
372907761f8SToby Isaac    Input Arguments:
373907761f8SToby Isaac +  q - the quadrature functional
374907761f8SToby Isaac .  imageDim - the dimension of the image of the transformation
375907761f8SToby Isaac .  origin - a point in the original space
376907761f8SToby Isaac .  originImage - the image of the origin under the transformation
377907761f8SToby Isaac .  J - the Jacobian of the image: an [imageDim x dim] matrix in row major order
37828222859SToby Isaac -  formDegree - transform the quadrature weights as k-forms of this form degree (if the number of components is a multiple of (dim choose formDegree), it is assumed that they represent multiple k-forms) [see PetscDTAltVPullback() for interpretation of formDegree]
379907761f8SToby Isaac 
380907761f8SToby Isaac    Output Arguments:
381907761f8SToby Isaac .  Jinvstarq - a quadrature rule where each point is the image of a point in the original quadrature rule, and where the k-form weights have been pulled-back by the pseudoinverse of J to the k-form weights in the image space.
382907761f8SToby Isaac 
383907761f8SToby Isaac    Note: the new quadrature rule will have a different number of components if spaces have different dimensions.  For example, pushing a 2-form forward from a two dimensional space to a three dimensional space changes the number of components from 1 to 3.
384907761f8SToby Isaac 
385907761f8SToby Isaac .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
386907761f8SToby Isaac @*/
38728222859SToby Isaac PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formDegree, PetscQuadrature *Jinvstarq)
388907761f8SToby Isaac {
389907761f8SToby Isaac   PetscInt         dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c;
390907761f8SToby Isaac   const PetscReal *points;
391907761f8SToby Isaac   const PetscReal *weights;
392907761f8SToby Isaac   PetscReal       *imagePoints, *imageWeights;
393907761f8SToby Isaac   PetscReal       *Jinv;
394907761f8SToby Isaac   PetscReal       *Jinvstar;
395907761f8SToby Isaac   PetscErrorCode   ierr;
396907761f8SToby Isaac 
397907761f8SToby Isaac   PetscFunctionBegin;
398907761f8SToby Isaac   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
39928222859SToby Isaac   if (imageDim < PetscAbsInt(formDegree)) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %D-form in %D dimensions", PetscAbsInt(formDegree), imageDim);
400907761f8SToby Isaac   ierr = PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights);CHKERRQ(ierr);
40128222859SToby Isaac   ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &formSize);CHKERRQ(ierr);
402907761f8SToby Isaac   if (Nc % formSize) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Number of components %D is not a multiple of formSize %D\n", Nc, formSize);
403907761f8SToby Isaac   Ncopies = Nc / formSize;
40428222859SToby Isaac   ierr = PetscDTBinomialInt(imageDim, PetscAbsInt(formDegree), &imageFormSize);CHKERRQ(ierr);
405907761f8SToby Isaac   imageNc = Ncopies * imageFormSize;
406907761f8SToby Isaac   ierr = PetscMalloc1(Npoints * imageDim, &imagePoints);CHKERRQ(ierr);
407907761f8SToby Isaac   ierr = PetscMalloc1(Npoints * imageNc, &imageWeights);CHKERRQ(ierr);
408907761f8SToby Isaac   ierr = PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar);CHKERRQ(ierr);
409907761f8SToby Isaac   ierr = PetscDTJacobianInverse_Internal(dim, imageDim, J, Jinv);CHKERRQ(ierr);
41028222859SToby Isaac   ierr = PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formDegree, Jinvstar);CHKERRQ(ierr);
411907761f8SToby Isaac   for (pt = 0; pt < Npoints; pt++) {
412907761f8SToby Isaac     const PetscReal *point = &points[pt * dim];
413907761f8SToby Isaac     PetscReal       *imagePoint = &imagePoints[pt * imageDim];
414907761f8SToby Isaac 
415907761f8SToby Isaac     for (i = 0; i < imageDim; i++) {
416907761f8SToby Isaac       PetscReal val = originImage[i];
417907761f8SToby Isaac 
418907761f8SToby Isaac       for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]);
419907761f8SToby Isaac       imagePoint[i] = val;
420907761f8SToby Isaac     }
421907761f8SToby Isaac     for (c = 0; c < Ncopies; c++) {
422907761f8SToby Isaac       const PetscReal *form = &weights[pt * Nc + c * formSize];
423907761f8SToby Isaac       PetscReal       *imageForm = &imageWeights[pt * imageNc + c * imageFormSize];
424907761f8SToby Isaac 
425907761f8SToby Isaac       for (i = 0; i < imageFormSize; i++) {
426907761f8SToby Isaac         PetscReal val = 0.;
427907761f8SToby Isaac 
428907761f8SToby Isaac         for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j];
429907761f8SToby Isaac         imageForm[i] = val;
430907761f8SToby Isaac       }
431907761f8SToby Isaac     }
432907761f8SToby Isaac   }
433907761f8SToby Isaac   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq);CHKERRQ(ierr);
434907761f8SToby Isaac   ierr = PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights);CHKERRQ(ierr);
435907761f8SToby Isaac   ierr = PetscFree2(Jinv, Jinvstar);CHKERRQ(ierr);
436907761f8SToby Isaac   PetscFunctionReturn(0);
437907761f8SToby Isaac }
438907761f8SToby Isaac 
43940d8ff71SMatthew G. Knepley /*@C
44040d8ff71SMatthew G. Knepley   PetscQuadratureSetData - Sets the data defining the quadrature
44140d8ff71SMatthew G. Knepley 
44240d8ff71SMatthew G. Knepley   Not collective
44340d8ff71SMatthew G. Knepley 
44440d8ff71SMatthew G. Knepley   Input Parameters:
44540d8ff71SMatthew G. Knepley + q  - The PetscQuadrature object
44640d8ff71SMatthew G. Knepley . dim - The spatial dimension
447e2b35d93SBarry Smith . Nc - The number of components
44840d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
44940d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
45040d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
45140d8ff71SMatthew G. Knepley 
452c99e0549SMatthew 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.
453f2fd9e53SMatthew G. Knepley 
45440d8ff71SMatthew G. Knepley   Level: intermediate
45540d8ff71SMatthew G. Knepley 
45640d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
45740d8ff71SMatthew G. Knepley @*/
458a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
45921454ff5SMatthew G. Knepley {
46021454ff5SMatthew G. Knepley   PetscFunctionBegin;
4612cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
46221454ff5SMatthew G. Knepley   if (dim >= 0)     q->dim       = dim;
463a6b92713SMatthew G. Knepley   if (Nc >= 0)      q->Nc        = Nc;
46421454ff5SMatthew G. Knepley   if (npoints >= 0) q->numPoints = npoints;
46521454ff5SMatthew G. Knepley   if (points) {
46621454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
46721454ff5SMatthew G. Knepley     q->points = points;
46821454ff5SMatthew G. Knepley   }
46921454ff5SMatthew G. Knepley   if (weights) {
47021454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
47121454ff5SMatthew G. Knepley     q->weights = weights;
47221454ff5SMatthew G. Knepley   }
473f9fd7fdbSMatthew G. Knepley   PetscFunctionReturn(0);
474f9fd7fdbSMatthew G. Knepley }
475f9fd7fdbSMatthew G. Knepley 
476d9bac1caSLisandro Dalcin static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v)
477d9bac1caSLisandro Dalcin {
478d9bac1caSLisandro Dalcin   PetscInt          q, d, c;
479d9bac1caSLisandro Dalcin   PetscViewerFormat format;
480d9bac1caSLisandro Dalcin   PetscErrorCode    ierr;
481d9bac1caSLisandro Dalcin 
482d9bac1caSLisandro Dalcin   PetscFunctionBegin;
483c74b4a09SMatthew G. Knepley   if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D) with %D components\n", quad->order, quad->numPoints, quad->dim, quad->Nc);CHKERRQ(ierr);}
484c74b4a09SMatthew G. Knepley   else              {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);}
485d9bac1caSLisandro Dalcin   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
486d9bac1caSLisandro Dalcin   if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
487d9bac1caSLisandro Dalcin   for (q = 0; q < quad->numPoints; ++q) {
488c74b4a09SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr);
489d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr);
490d9bac1caSLisandro Dalcin     for (d = 0; d < quad->dim; ++d) {
491d9bac1caSLisandro Dalcin       if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
492d9bac1caSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
493d9bac1caSLisandro Dalcin     }
494d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr);
495c74b4a09SMatthew G. Knepley     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);}
496d9bac1caSLisandro Dalcin     for (c = 0; c < quad->Nc; ++c) {
497d9bac1caSLisandro Dalcin       if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
498c74b4a09SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr);
499d9bac1caSLisandro Dalcin     }
500d9bac1caSLisandro Dalcin     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);}
501d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr);
502d9bac1caSLisandro Dalcin     ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr);
503d9bac1caSLisandro Dalcin   }
504d9bac1caSLisandro Dalcin   PetscFunctionReturn(0);
505d9bac1caSLisandro Dalcin }
506d9bac1caSLisandro Dalcin 
50740d8ff71SMatthew G. Knepley /*@C
50840d8ff71SMatthew G. Knepley   PetscQuadratureView - Views a PetscQuadrature object
50940d8ff71SMatthew G. Knepley 
510d083f849SBarry Smith   Collective on quad
51140d8ff71SMatthew G. Knepley 
51240d8ff71SMatthew G. Knepley   Input Parameters:
513d9bac1caSLisandro Dalcin + quad  - The PetscQuadrature object
51440d8ff71SMatthew G. Knepley - viewer - The PetscViewer object
51540d8ff71SMatthew G. Knepley 
51640d8ff71SMatthew G. Knepley   Level: beginner
51740d8ff71SMatthew G. Knepley 
51840d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
51940d8ff71SMatthew G. Knepley @*/
520f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
521f9fd7fdbSMatthew G. Knepley {
522d9bac1caSLisandro Dalcin   PetscBool      iascii;
523f9fd7fdbSMatthew G. Knepley   PetscErrorCode ierr;
524f9fd7fdbSMatthew G. Knepley 
525f9fd7fdbSMatthew G. Knepley   PetscFunctionBegin;
526d9bac1caSLisandro Dalcin   PetscValidHeader(quad, 1);
527d9bac1caSLisandro Dalcin   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
528d9bac1caSLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);}
529d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
530d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
531d9bac1caSLisandro Dalcin   if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);}
532d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
533bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
534bfa639d9SMatthew G. Knepley }
535bfa639d9SMatthew G. Knepley 
53689710940SMatthew G. Knepley /*@C
53789710940SMatthew G. Knepley   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
53889710940SMatthew G. Knepley 
53989710940SMatthew G. Knepley   Not collective
54089710940SMatthew G. Knepley 
54189710940SMatthew G. Knepley   Input Parameter:
54289710940SMatthew G. Knepley + q - The original PetscQuadrature
54389710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into
54489710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement
54589710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement
54689710940SMatthew G. Knepley 
54789710940SMatthew G. Knepley   Output Parameters:
54889710940SMatthew G. Knepley . dim - The dimension
54989710940SMatthew G. Knepley 
55089710940SMatthew G. Knepley   Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
55189710940SMatthew G. Knepley 
552f5f57ec0SBarry Smith  Not available from Fortran
553f5f57ec0SBarry Smith 
55489710940SMatthew G. Knepley   Level: intermediate
55589710940SMatthew G. Knepley 
55689710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
55789710940SMatthew G. Knepley @*/
55889710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
55989710940SMatthew G. Knepley {
56089710940SMatthew G. Knepley   const PetscReal *points,    *weights;
56189710940SMatthew G. Knepley   PetscReal       *pointsRef, *weightsRef;
562a6b92713SMatthew G. Knepley   PetscInt         dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
56389710940SMatthew G. Knepley   PetscErrorCode   ierr;
56489710940SMatthew G. Knepley 
56589710940SMatthew G. Knepley   PetscFunctionBegin;
5662cd22861SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
56789710940SMatthew G. Knepley   PetscValidPointer(v0, 3);
56889710940SMatthew G. Knepley   PetscValidPointer(jac, 4);
56989710940SMatthew G. Knepley   PetscValidPointer(qref, 5);
57089710940SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr);
57189710940SMatthew G. Knepley   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
572a6b92713SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr);
57389710940SMatthew G. Knepley   npointsRef = npoints*numSubelements;
57489710940SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr);
575a6b92713SMatthew G. Knepley   ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr);
57689710940SMatthew G. Knepley   for (c = 0; c < numSubelements; ++c) {
57789710940SMatthew G. Knepley     for (p = 0; p < npoints; ++p) {
57889710940SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
57989710940SMatthew G. Knepley         pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
58089710940SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
58189710940SMatthew G. Knepley           pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
58289710940SMatthew G. Knepley         }
58389710940SMatthew G. Knepley       }
58489710940SMatthew G. Knepley       /* Could also use detJ here */
585a6b92713SMatthew G. Knepley       for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements;
58689710940SMatthew G. Knepley     }
58789710940SMatthew G. Knepley   }
58889710940SMatthew G. Knepley   ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr);
589a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr);
59089710940SMatthew G. Knepley   PetscFunctionReturn(0);
59189710940SMatthew G. Knepley }
59289710940SMatthew G. Knepley 
59337045ce4SJed Brown /*@
59437045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
59537045ce4SJed Brown 
59637045ce4SJed Brown    Not Collective
59737045ce4SJed Brown 
59837045ce4SJed Brown    Input Arguments:
59937045ce4SJed Brown +  npoints - number of spatial points to evaluate at
60037045ce4SJed Brown .  points - array of locations to evaluate at
60137045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
60237045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
60337045ce4SJed Brown 
60437045ce4SJed Brown    Output Arguments:
6050298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
6060298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
6070298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
60837045ce4SJed Brown 
60937045ce4SJed Brown    Level: intermediate
61037045ce4SJed Brown 
61137045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
61237045ce4SJed Brown @*/
61337045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
61437045ce4SJed Brown {
61537045ce4SJed Brown   PetscInt i,maxdegree;
61637045ce4SJed Brown 
61737045ce4SJed Brown   PetscFunctionBegin;
61837045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
61937045ce4SJed Brown   maxdegree = degrees[ndegree-1];
62037045ce4SJed Brown   for (i=0; i<npoints; i++) {
62137045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
62237045ce4SJed Brown     PetscInt  j,k;
62337045ce4SJed Brown     x    = points[i];
62437045ce4SJed Brown     pm2  = 0;
62537045ce4SJed Brown     pm1  = 1;
62637045ce4SJed Brown     pd2  = 0;
62737045ce4SJed Brown     pd1  = 0;
62837045ce4SJed Brown     pdd2 = 0;
62937045ce4SJed Brown     pdd1 = 0;
63037045ce4SJed Brown     k    = 0;
63137045ce4SJed Brown     if (degrees[k] == 0) {
63237045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
63337045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
63437045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
63537045ce4SJed Brown       k++;
63637045ce4SJed Brown     }
63737045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
63837045ce4SJed Brown       PetscReal p,d,dd;
63937045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
64037045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
64137045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
64237045ce4SJed Brown       pm2  = pm1;
64337045ce4SJed Brown       pm1  = p;
64437045ce4SJed Brown       pd2  = pd1;
64537045ce4SJed Brown       pd1  = d;
64637045ce4SJed Brown       pdd2 = pdd1;
64737045ce4SJed Brown       pdd1 = dd;
64837045ce4SJed Brown       if (degrees[k] == j) {
64937045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
65037045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
65137045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
65237045ce4SJed Brown       }
65337045ce4SJed Brown     }
65437045ce4SJed Brown   }
65537045ce4SJed Brown   PetscFunctionReturn(0);
65637045ce4SJed Brown }
65737045ce4SJed Brown 
658*e6a796c3SToby Isaac /* if we can cobble together an eigenvalue / eigenvector routine for a symmetric tridiagonal system, we should use
659*e6a796c3SToby Isaac  * Golub & Welsch (Gauss-Jacobi) or Golub (Gauss-Lobatto-Jacobi) to compute Gaussian quadrature */
660*e6a796c3SToby Isaac #if (!defined(PETSC_MISSING_LAPACK_STEQR) || !defined(PETSC_MISSING_LAPACK_STEGR))
661*e6a796c3SToby Isaac #define PETSCDTGAUSSIANQUADRATURE_EIG 1
662*e6a796c3SToby Isaac #endif
663*e6a796c3SToby Isaac 
664*e6a796c3SToby Isaac /* solve the symmetric tridiagonal eigenvalue system, writing the eigenvalues into eigs and the eigenvectors into V
665*e6a796c3SToby Isaac  * with lds n; diag and subdiag are overwritten */
666*e6a796c3SToby Isaac static PetscErrorCode PetscDTSymmetricTridiagonalEigensolve(PetscInt n, PetscReal diag[], PetscReal subdiag[],
667*e6a796c3SToby Isaac                                                             PetscReal eigs[], PetscScalar V[])
668*e6a796c3SToby Isaac {
669*e6a796c3SToby Isaac   char jobz = 'V'; /* eigenvalues and eigenvectors */
670*e6a796c3SToby Isaac   char range = 'A'; /* all eigenvalues will be found */
671*e6a796c3SToby Isaac   PetscReal VL = 0.; /* ignored because range is 'A' */
672*e6a796c3SToby Isaac   PetscReal VU = 0.; /* ignored because range is 'A' */
673*e6a796c3SToby Isaac   PetscBLASInt IL = 0; /* ignored because range is 'A' */
674*e6a796c3SToby Isaac   PetscBLASInt IU = 0; /* ignored because range is 'A' */
675*e6a796c3SToby Isaac   PetscReal abstol = 0.; /* unused */
676*e6a796c3SToby Isaac   PetscBLASInt bn, bm, ldz; /* bm will equal bn on exit */
677*e6a796c3SToby Isaac   PetscBLASInt *isuppz;
678*e6a796c3SToby Isaac   PetscBLASInt lwork, liwork;
679*e6a796c3SToby Isaac   PetscReal workquery;
680*e6a796c3SToby Isaac   PetscBLASInt  iworkquery;
681*e6a796c3SToby Isaac   PetscBLASInt *iwork;
682*e6a796c3SToby Isaac   PetscBLASInt info;
683*e6a796c3SToby Isaac   PetscReal *work = NULL;
684*e6a796c3SToby Isaac   PetscErrorCode ierr;
685*e6a796c3SToby Isaac 
686*e6a796c3SToby Isaac   PetscFunctionBegin;
687*e6a796c3SToby Isaac #if !defined(PETSCDTGAUSSIANQUADRATURE_EIG)
688*e6a796c3SToby Isaac   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
689*e6a796c3SToby Isaac #endif
690*e6a796c3SToby Isaac   ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr);
691*e6a796c3SToby Isaac   ierr = PetscBLASIntCast(n, &ldz);CHKERRQ(ierr);
692*e6a796c3SToby Isaac #if !defined(PETSC_MISSING_LAPACK_STEGR)
693*e6a796c3SToby Isaac   ierr = PetscMalloc1(2 * n, &isuppz);CHKERRQ(ierr);
694*e6a796c3SToby Isaac   lwork = -1;
695*e6a796c3SToby Isaac   liwork = -1;
696*e6a796c3SToby Isaac   PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,&workquery,&lwork,&iworkquery,&liwork,&info));
697*e6a796c3SToby Isaac   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error");
698*e6a796c3SToby Isaac   lwork = (PetscBLASInt) workquery;
699*e6a796c3SToby Isaac   liwork = (PetscBLASInt) iworkquery;
700*e6a796c3SToby Isaac   ierr = PetscMalloc2(lwork, &work, liwork, &iwork);CHKERRQ(ierr);
701*e6a796c3SToby Isaac   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
702*e6a796c3SToby Isaac   PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,work,&lwork,iwork,&liwork,&info));
703*e6a796c3SToby Isaac   ierr = PetscFPTrapPop();CHKERRQ(ierr);
704*e6a796c3SToby Isaac   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error");
705*e6a796c3SToby Isaac   ierr = PetscFree2(work, iwork);CHKERRQ(ierr);
706*e6a796c3SToby Isaac   ierr = PetscFree(isuppz);CHKERRQ(ierr);
707*e6a796c3SToby Isaac #elif !defined(PETSC_MISSING_LAPACK_STEQR)
708*e6a796c3SToby Isaac   jobz = 'I'; /* Compute eigenvalues and eigenvectors of the
709*e6a796c3SToby Isaac                  tridiagonal matrix.  Z is initialized to the identity
710*e6a796c3SToby Isaac                  matrix. */
711*e6a796c3SToby Isaac   ierr = PetscMalloc1(PetscMax(1,2*n-2),&work);CHKERRQ(ierr);
712*e6a796c3SToby Isaac   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&bn,diag,subdiag,V,&ldz,work,&info));
713*e6a796c3SToby Isaac   ierr = PetscFPTrapPop();CHKERRQ(ierr);
714*e6a796c3SToby Isaac   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
715*e6a796c3SToby Isaac   ierr = PetscFree(work);CHKERRQ(ierr);
716*e6a796c3SToby Isaac   ierr = PetscArraycpy(eigs,diag,n);CHKERRQ(ierr);
717*e6a796c3SToby Isaac #endif
718*e6a796c3SToby Isaac   PetscFunctionReturn(0);
719*e6a796c3SToby Isaac }
720*e6a796c3SToby Isaac 
721*e6a796c3SToby Isaac /* Formula for the weights at the endpoints (-1 and 1) of Gauss-Lobatto-Jacobi
722*e6a796c3SToby Isaac  * quadrature rules on the interval [-1, 1] */
723*e6a796c3SToby Isaac static PetscErrorCode PetscDTGaussLobattoJacobiEndweights_Internal(PetscInt n, PetscReal alpha, PetscReal beta, PetscReal *leftw, PetscReal *rightw)
724*e6a796c3SToby Isaac {
725*e6a796c3SToby Isaac   PetscReal twoab1;
726*e6a796c3SToby Isaac   PetscInt  m = n - 2;
727*e6a796c3SToby Isaac   PetscReal a = alpha + 1.;
728*e6a796c3SToby Isaac   PetscReal b = beta + 1.;
729*e6a796c3SToby Isaac   PetscReal gra, grb;
730*e6a796c3SToby Isaac   PetscErrorCode ierr;
731*e6a796c3SToby Isaac 
732*e6a796c3SToby Isaac   PetscFunctionBegin;
733*e6a796c3SToby Isaac   twoab1 = PetscPowReal(2., a + b - 1.);
734*e6a796c3SToby Isaac #if defined(PETSC_HAVE_LGAMMA)
735*e6a796c3SToby Isaac   grb = PetscExpReal(2. * PetscLGamma(b+1.) + PetscLGamma(m+1.) + PetscLGamma(m+a+1.) -
736*e6a796c3SToby Isaac                      (PetscLGamma(m+b+1) + PetscLGamma(m+a+b+1.)));
737*e6a796c3SToby Isaac   gra = PetscExpReal(2. * PetscLGamma(a+1.) + PetscLGamma(m+1.) + PetscLGamma(m+b+1.) -
738*e6a796c3SToby Isaac                      (PetscLGamma(m+a+1) + PetscLGamma(m+a+b+1.)));
739*e6a796c3SToby Isaac #else
740*e6a796c3SToby Isaac   {
741*e6a796c3SToby Isaac     PetscInt alphai = (PetscInt) alpha;
742*e6a796c3SToby Isaac     PetscInt betai = (PetscInt) beta;
743*e6a796c3SToby Isaac 
744*e6a796c3SToby Isaac     if ((PetscReal) alphai == alpha && (PetscReal) betai == beta) {
745*e6a796c3SToby Isaac       PetscReal binom1, binom2;
746*e6a796c3SToby Isaac 
747*e6a796c3SToby Isaac       ierr = PetscDTBinomial(m+b, b, &binom1);CHKERRQ(ierr);
748*e6a796c3SToby Isaac       ierr = PetscDTBinomial(m+a+b, b, &binom2);CHKERRQ(ierr);
749*e6a796c3SToby Isaac       grb = 1./ (binom1 * binom2);
750*e6a796c3SToby Isaac       ierr = PetscDTBinomial(m+a, a, &binom1);CHKERRQ(ierr);
751*e6a796c3SToby Isaac       ierr = PetscDTBinomial(m+a+b, a, &binom2);CHKERRQ(ierr);
752*e6a796c3SToby Isaac       gra = 1./ (binom1 * binom2);
753*e6a796c3SToby Isaac     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable.");
754*e6a796c3SToby Isaac   }
755*e6a796c3SToby Isaac #endif
756*e6a796c3SToby Isaac   *leftw = twoab1 * grb / b;
757*e6a796c3SToby Isaac   *rightw = twoab1 * gra / a;
758*e6a796c3SToby Isaac   PetscFunctionReturn(0);
759*e6a796c3SToby Isaac }
760*e6a796c3SToby Isaac 
761*e6a796c3SToby Isaac /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
762*e6a796c3SToby Isaac    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
763*e6a796c3SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
764*e6a796c3SToby Isaac {
765*e6a796c3SToby Isaac   PetscReal apb, pn1, pn2;
766*e6a796c3SToby Isaac   PetscInt  k;
767*e6a796c3SToby Isaac 
768*e6a796c3SToby Isaac   PetscFunctionBegin;
769*e6a796c3SToby Isaac   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
770*e6a796c3SToby Isaac   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
771*e6a796c3SToby Isaac   apb = a + b;
772*e6a796c3SToby Isaac   pn2 = 1.0;
773*e6a796c3SToby Isaac   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
774*e6a796c3SToby Isaac   *P  = 0.0;
775*e6a796c3SToby Isaac   for (k = 2; k < n+1; ++k) {
776*e6a796c3SToby Isaac     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
777*e6a796c3SToby Isaac     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
778*e6a796c3SToby Isaac     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
779*e6a796c3SToby Isaac     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
780*e6a796c3SToby Isaac 
781*e6a796c3SToby Isaac     a2  = a2 / a1;
782*e6a796c3SToby Isaac     a3  = a3 / a1;
783*e6a796c3SToby Isaac     a4  = a4 / a1;
784*e6a796c3SToby Isaac     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
785*e6a796c3SToby Isaac     pn2 = pn1;
786*e6a796c3SToby Isaac     pn1 = *P;
787*e6a796c3SToby Isaac   }
788*e6a796c3SToby Isaac   PetscFunctionReturn(0);
789*e6a796c3SToby Isaac }
790*e6a796c3SToby Isaac 
791*e6a796c3SToby Isaac /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
792*e6a796c3SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscInt k, PetscReal *P)
793*e6a796c3SToby Isaac {
794*e6a796c3SToby Isaac   PetscReal      nP;
795*e6a796c3SToby Isaac   PetscInt       i;
796*e6a796c3SToby Isaac   PetscErrorCode ierr;
797*e6a796c3SToby Isaac 
798*e6a796c3SToby Isaac   PetscFunctionBegin;
799*e6a796c3SToby Isaac   if (k > n) {*P = 0.0; PetscFunctionReturn(0);}
800*e6a796c3SToby Isaac   ierr = PetscDTComputeJacobi(a+k, b+k, n-k, x, &nP);CHKERRQ(ierr);
801*e6a796c3SToby Isaac   for (i = 0; i < k; i++) nP *= (a + b + n + 1. + i) * 0.5;
802*e6a796c3SToby Isaac   *P = nP;
803*e6a796c3SToby Isaac   PetscFunctionReturn(0);
804*e6a796c3SToby Isaac }
805*e6a796c3SToby Isaac 
806*e6a796c3SToby Isaac /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
807*e6a796c3SToby Isaac PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
808*e6a796c3SToby Isaac {
809*e6a796c3SToby Isaac   PetscFunctionBegin;
810*e6a796c3SToby Isaac   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
811*e6a796c3SToby Isaac   *eta = y;
812*e6a796c3SToby Isaac   PetscFunctionReturn(0);
813*e6a796c3SToby Isaac }
814*e6a796c3SToby Isaac 
815*e6a796c3SToby Isaac static PetscErrorCode PetscDTGaussJacobiQuadrature_Newton_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[])
816*e6a796c3SToby Isaac {
817*e6a796c3SToby Isaac   PetscInt       maxIter = 100;
818*e6a796c3SToby Isaac   PetscReal      eps     = 1.0e-8;
819*e6a796c3SToby Isaac   PetscReal      a1, a2, a3, a4, a5, a6;
820*e6a796c3SToby Isaac   PetscInt       k;
821*e6a796c3SToby Isaac   PetscErrorCode ierr;
822*e6a796c3SToby Isaac 
823*e6a796c3SToby Isaac   PetscFunctionBegin;
824*e6a796c3SToby Isaac 
825*e6a796c3SToby Isaac   a1      = PetscPowReal(2.0, a+b+1);
826*e6a796c3SToby Isaac #if defined(PETSC_HAVE_TGAMMA)
827*e6a796c3SToby Isaac   a2      = PetscTGamma(a + npoints + 1);
828*e6a796c3SToby Isaac   a3      = PetscTGamma(b + npoints + 1);
829*e6a796c3SToby Isaac   a4      = PetscTGamma(a + b + npoints + 1);
830*e6a796c3SToby Isaac #else
831*e6a796c3SToby Isaac   {
832*e6a796c3SToby Isaac     PetscInt ia, ib;
833*e6a796c3SToby Isaac 
834*e6a796c3SToby Isaac     ia = (PetscInt) a;
835*e6a796c3SToby Isaac     ib = (PetscInt) b;
836*e6a796c3SToby 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 */
837*e6a796c3SToby Isaac       ierr = PetscDTFactorial(ia + npoints, &a2);CHKERRQ(ierr);
838*e6a796c3SToby Isaac       ierr = PetscDTFactorial(ib + npoints, &a3);CHKERRQ(ierr);
839*e6a796c3SToby Isaac       ierr = PetscDTFactorial(ia + ib + npoints, &a4);CHKERRQ(ierr);
840*e6a796c3SToby Isaac     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
841*e6a796c3SToby Isaac   }
842*e6a796c3SToby Isaac #endif
843*e6a796c3SToby Isaac 
844*e6a796c3SToby Isaac   ierr = PetscDTFactorial(npoints, &a5);CHKERRQ(ierr);
845*e6a796c3SToby Isaac   a6   = a1 * a2 * a3 / a4 / a5;
846*e6a796c3SToby Isaac   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
847*e6a796c3SToby Isaac    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
848*e6a796c3SToby Isaac   for (k = 0; k < npoints; ++k) {
849*e6a796c3SToby Isaac     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
850*e6a796c3SToby Isaac     PetscInt  j;
851*e6a796c3SToby Isaac 
852*e6a796c3SToby Isaac     if (k > 0) r = 0.5 * (r + x[k-1]);
853*e6a796c3SToby Isaac     for (j = 0; j < maxIter; ++j) {
854*e6a796c3SToby Isaac       PetscReal s = 0.0, delta, f, fp;
855*e6a796c3SToby Isaac       PetscInt  i;
856*e6a796c3SToby Isaac 
857*e6a796c3SToby Isaac       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
858*e6a796c3SToby Isaac       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
859*e6a796c3SToby Isaac       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, 1, &fp);CHKERRQ(ierr);
860*e6a796c3SToby Isaac       delta = f / (fp - f * s);
861*e6a796c3SToby Isaac       r     = r - delta;
862*e6a796c3SToby Isaac       if (PetscAbsReal(delta) < eps) break;
863*e6a796c3SToby Isaac     }
864*e6a796c3SToby Isaac     x[k] = r;
865*e6a796c3SToby Isaac     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], 1, &dP);CHKERRQ(ierr);
866*e6a796c3SToby Isaac     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
867*e6a796c3SToby Isaac   }
868*e6a796c3SToby Isaac   PetscFunctionReturn(0);
869*e6a796c3SToby Isaac }
870*e6a796c3SToby Isaac 
871*e6a796c3SToby Isaac /* Compute the diagonals of the Jacobi matrix used in Golub (& Welsch) algorithms for Gauss-(Lobatto)-Jacobi
872*e6a796c3SToby Isaac  * quadrature weight calculations on [-1,1] for exponents (1. + x)^a (1.-x)^b */
873*e6a796c3SToby Isaac static PetscErrorCode PetscDTJacobiMatrix_Internal(PetscInt nPoints, PetscReal a, PetscReal b, PetscReal *d, PetscReal *s)
874*e6a796c3SToby Isaac {
875*e6a796c3SToby Isaac   PetscInt       i;
876*e6a796c3SToby Isaac 
877*e6a796c3SToby Isaac   PetscFunctionBegin;
878*e6a796c3SToby Isaac   for (i = 0; i < nPoints; i++) {
879*e6a796c3SToby Isaac     PetscReal A, B, Ap, C;
880*e6a796c3SToby Isaac     PetscInt  si = i+1;
881*e6a796c3SToby Isaac 
882*e6a796c3SToby Isaac     /* Jacobi three term recurrence */
883*e6a796c3SToby Isaac     if (si > 1) {
884*e6a796c3SToby Isaac       A = (2.*si+a+b)*(2*si+a+b-1.) / (2.*si*(si+a+b));
885*e6a796c3SToby Isaac       B = (a*a-b*b)*(2*si+a+b-1.) / (2.*si*(si+a+b)*(2*si+a+b-2));
886*e6a796c3SToby Isaac     } else {
887*e6a796c3SToby Isaac       A = (a+b+2.) / 2.;
888*e6a796c3SToby Isaac       B = (a-b) / 2.;
889*e6a796c3SToby Isaac     }
890*e6a796c3SToby Isaac     Ap = (2*(si+1)+a+b)*(2*(si+1)+a+b-1) / (2*(si+1)*(si+1+a+b));
891*e6a796c3SToby Isaac     C = (2*((si+1)+a-1)*((si+1)+b-1)*(2*(si+1)+a+b) / (2*(si+1)*((si+1)+a+b)*(2*(si+1)+a+b-2)));
892*e6a796c3SToby Isaac     d[i] = -B / A;
893*e6a796c3SToby Isaac     s[i] = (C / (A * Ap));
894*e6a796c3SToby Isaac   }
895*e6a796c3SToby Isaac   PetscFunctionReturn(0);
896*e6a796c3SToby Isaac }
897*e6a796c3SToby Isaac 
898*e6a796c3SToby Isaac static PetscErrorCode PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
899*e6a796c3SToby Isaac {
900*e6a796c3SToby Isaac   PetscReal mu0;
901*e6a796c3SToby Isaac   PetscReal ga, gb, gab;
902*e6a796c3SToby Isaac   PetscInt i;
903*e6a796c3SToby Isaac   PetscErrorCode ierr;
904*e6a796c3SToby Isaac 
905*e6a796c3SToby Isaac   PetscFunctionBegin;
906*e6a796c3SToby Isaac   ierr = PetscCitationsRegister(GolubWelschCitation, &GolubWelschCite);CHKERRQ(ierr);
907*e6a796c3SToby Isaac 
908*e6a796c3SToby Isaac #if defined(PETSC_HAVE_TGAMMA)
909*e6a796c3SToby Isaac   ga  = PetscTGamma(a + 1);
910*e6a796c3SToby Isaac   gb  = PetscTGamma(b + 1);
911*e6a796c3SToby Isaac   gab = PetscTGamma(a + b + 2);
912*e6a796c3SToby Isaac #else
913*e6a796c3SToby Isaac   {
914*e6a796c3SToby Isaac     PetscInt ia, ib;
915*e6a796c3SToby Isaac 
916*e6a796c3SToby Isaac     ia = (PetscInt) a;
917*e6a796c3SToby Isaac     ib = (PetscInt) b;
918*e6a796c3SToby Isaac     if (ia == a && ib == b && ia + 1 > 0 && ib + 1 > 0 && ia + ib + 2 > 0) { /* All gamma(x) terms are (x-1)! terms */
919*e6a796c3SToby Isaac       ierr = PetscDTFactorial(ia, &ga);CHKERRQ(ierr);
920*e6a796c3SToby Isaac       ierr = PetscDTFactorial(ib, &gb);CHKERRQ(ierr);
921*e6a796c3SToby Isaac       ierr = PetscDTFactorial(ia + ib + 1, &gb);CHKERRQ(ierr);
922*e6a796c3SToby Isaac     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
923*e6a796c3SToby Isaac   }
924*e6a796c3SToby Isaac #endif
925*e6a796c3SToby Isaac   mu0 = PetscPowReal(2.,a + b + 1.) * ga * gb / gab;
926*e6a796c3SToby Isaac 
927*e6a796c3SToby Isaac #if defined(PETSCDTGAUSSIANQUADRATURE_EIG)
928*e6a796c3SToby Isaac   {
929*e6a796c3SToby Isaac     PetscReal *diag, *subdiag;
930*e6a796c3SToby Isaac     PetscScalar *V;
931*e6a796c3SToby Isaac 
932*e6a796c3SToby Isaac     ierr = PetscMalloc2(npoints, &diag, npoints, &subdiag);CHKERRQ(ierr);
933*e6a796c3SToby Isaac     ierr = PetscMalloc1(npoints*npoints, &V);CHKERRQ(ierr);
934*e6a796c3SToby Isaac     ierr = PetscDTJacobiMatrix_Internal(npoints, a, b, diag, subdiag);CHKERRQ(ierr);
935*e6a796c3SToby Isaac     for (i = 0; i < npoints - 1; i++) subdiag[i] = PetscSqrtReal(subdiag[i]);
936*e6a796c3SToby Isaac     ierr = PetscDTSymmetricTridiagonalEigensolve(npoints, diag, subdiag, x, V);CHKERRQ(ierr);
937*e6a796c3SToby Isaac     for (i = 0; i < npoints; i++) w[i] = PetscSqr(V[i * npoints]) * mu0;
938*e6a796c3SToby Isaac     ierr = PetscFree(V);CHKERRQ(ierr);
939*e6a796c3SToby Isaac     ierr = PetscFree2(diag, subdiag);CHKERRQ(ierr);
940*e6a796c3SToby Isaac   }
941*e6a796c3SToby Isaac #else
942*e6a796c3SToby Isaac   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
943*e6a796c3SToby Isaac #endif
944*e6a796c3SToby Isaac   PetscFunctionReturn(0);
945*e6a796c3SToby Isaac }
946*e6a796c3SToby Isaac 
947*e6a796c3SToby Isaac static PetscErrorCode PetscDTGaussJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
948*e6a796c3SToby Isaac {
949*e6a796c3SToby Isaac   PetscErrorCode ierr;
950*e6a796c3SToby Isaac 
951*e6a796c3SToby Isaac   PetscFunctionBegin;
952*e6a796c3SToby Isaac   if (npoints < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive");
953*e6a796c3SToby Isaac   /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
954*e6a796c3SToby Isaac   if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1.");
955*e6a796c3SToby Isaac   if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1.");
956*e6a796c3SToby Isaac 
957*e6a796c3SToby Isaac   if (newton) {
958*e6a796c3SToby Isaac     ierr = PetscDTGaussJacobiQuadrature_Newton_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr);
959*e6a796c3SToby Isaac   } else {
960*e6a796c3SToby Isaac     ierr = PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr);
961*e6a796c3SToby Isaac   }
962*e6a796c3SToby Isaac   if (alpha == beta) { /* symmetrize */
963*e6a796c3SToby Isaac     PetscInt i;
964*e6a796c3SToby Isaac     for (i = 0; i < (npoints + 1) / 2; i++) {
965*e6a796c3SToby Isaac       PetscInt  j  = npoints - 1 - i;
966*e6a796c3SToby Isaac       PetscReal xi = x[i];
967*e6a796c3SToby Isaac       PetscReal xj = x[j];
968*e6a796c3SToby Isaac       PetscReal wi = w[i];
969*e6a796c3SToby Isaac       PetscReal wj = w[j];
970*e6a796c3SToby Isaac 
971*e6a796c3SToby Isaac       x[i] = (xi - xj) / 2.;
972*e6a796c3SToby Isaac       x[j] = (xj - xi) / 2.;
973*e6a796c3SToby Isaac       w[i] = w[j] = (wi + wj) / 2.;
974*e6a796c3SToby Isaac     }
975*e6a796c3SToby Isaac   }
976*e6a796c3SToby Isaac   PetscFunctionReturn(0);
977*e6a796c3SToby Isaac }
978*e6a796c3SToby Isaac 
979*e6a796c3SToby Isaac PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
980*e6a796c3SToby Isaac {
981*e6a796c3SToby Isaac   PetscErrorCode ierr;
982*e6a796c3SToby Isaac 
983*e6a796c3SToby Isaac   PetscFunctionBegin;
984*e6a796c3SToby Isaac   ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, alpha, beta, x, w, gaussQuadratureNewton);CHKERRQ(ierr);
985*e6a796c3SToby Isaac   PetscFunctionReturn(0);
986*e6a796c3SToby Isaac }
987*e6a796c3SToby Isaac 
988*e6a796c3SToby Isaac static PetscErrorCode PetscDTGaussLobattoJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
989*e6a796c3SToby Isaac {
990*e6a796c3SToby Isaac   PetscInt       i;
991*e6a796c3SToby Isaac   PetscErrorCode ierr;
992*e6a796c3SToby Isaac 
993*e6a796c3SToby Isaac   PetscFunctionBegin;
994*e6a796c3SToby Isaac   if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive");
995*e6a796c3SToby Isaac   /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
996*e6a796c3SToby Isaac   if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1.");
997*e6a796c3SToby Isaac   if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1.");
998*e6a796c3SToby Isaac 
999*e6a796c3SToby Isaac   x[0] = -1.;
1000*e6a796c3SToby Isaac   x[npoints-1] = 1.;
1001*e6a796c3SToby Isaac   ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, alpha+1., beta+1., &x[1], &w[1], newton);CHKERRQ(ierr);
1002*e6a796c3SToby Isaac   for (i = 1; i < npoints - 1; i++) {
1003*e6a796c3SToby Isaac     w[i] /= (1. - x[i]*x[i]);
1004*e6a796c3SToby Isaac   }
1005*e6a796c3SToby Isaac   ierr = PetscDTGaussLobattoJacobiEndweights_Internal(npoints, alpha, beta, &w[0], &w[npoints-1]);CHKERRQ(ierr);
1006*e6a796c3SToby Isaac   PetscFunctionReturn(0);
1007*e6a796c3SToby Isaac }
1008*e6a796c3SToby Isaac 
100937045ce4SJed Brown /*@
1010*e6a796c3SToby Isaac    PetscDTGaussQuadrature - create Gauss-Legendre quadrature
101137045ce4SJed Brown 
101237045ce4SJed Brown    Not Collective
101337045ce4SJed Brown 
101437045ce4SJed Brown    Input Arguments:
101537045ce4SJed Brown +  npoints - number of points
101637045ce4SJed Brown .  a - left end of interval (often-1)
101737045ce4SJed Brown -  b - right end of interval (often +1)
101837045ce4SJed Brown 
101937045ce4SJed Brown    Output Arguments:
102037045ce4SJed Brown +  x - quadrature points
102137045ce4SJed Brown -  w - quadrature weights
102237045ce4SJed Brown 
102337045ce4SJed Brown    Level: intermediate
102437045ce4SJed Brown 
102537045ce4SJed Brown    References:
102696a0c994SBarry Smith .   1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
102737045ce4SJed Brown 
102837045ce4SJed Brown .seealso: PetscDTLegendreEval()
102937045ce4SJed Brown @*/
103037045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
103137045ce4SJed Brown {
103237045ce4SJed Brown   PetscInt       i;
1033*e6a796c3SToby Isaac   PetscErrorCode ierr;
103437045ce4SJed Brown 
103537045ce4SJed Brown   PetscFunctionBegin;
1036*e6a796c3SToby Isaac   ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, 0., 0., x, w, gaussQuadratureNewton);CHKERRQ(ierr);
1037*e6a796c3SToby Isaac   if (a != -1. || b != 1.) {
103837045ce4SJed Brown     for (i = 0; i < npoints; i++) {
1039*e6a796c3SToby Isaac       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1040*e6a796c3SToby Isaac       w[i] *= (b - a) / 2.;
104137045ce4SJed Brown     }
104237045ce4SJed Brown   }
104337045ce4SJed Brown   PetscFunctionReturn(0);
104437045ce4SJed Brown }
1045194825f6SJed Brown 
10468272889dSSatish Balay /*@C
10478272889dSSatish Balay    PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre
10488272889dSSatish Balay                       nodes of a given size on the domain [-1,1]
10498272889dSSatish Balay 
10508272889dSSatish Balay    Not Collective
10518272889dSSatish Balay 
10528272889dSSatish Balay    Input Parameter:
10538272889dSSatish Balay +  n - number of grid nodes
1054f2e8fe4dShannah_mairs -  type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
10558272889dSSatish Balay 
10568272889dSSatish Balay    Output Arguments:
10578272889dSSatish Balay +  x - quadrature points
10588272889dSSatish Balay -  w - quadrature weights
10598272889dSSatish Balay 
10608272889dSSatish Balay    Notes:
10618272889dSSatish Balay     For n > 30  the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
10628272889dSSatish Balay           close enough to the desired solution
10638272889dSSatish Balay 
10648272889dSSatish Balay    These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes
10658272889dSSatish Balay 
1066a8d69d7bSBarry Smith    See  https://epubs.siam.org/doi/abs/10.1137/110855442  https://epubs.siam.org/doi/abs/10.1137/120889873 for better ways to compute GLL nodes
10678272889dSSatish Balay 
10688272889dSSatish Balay    Level: intermediate
10698272889dSSatish Balay 
10708272889dSSatish Balay .seealso: PetscDTGaussQuadrature()
10718272889dSSatish Balay 
10728272889dSSatish Balay @*/
1073916e780bShannah_mairs PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w)
10748272889dSSatish Balay {
1075*e6a796c3SToby Isaac   PetscBool      newton;
10768272889dSSatish Balay   PetscErrorCode ierr;
10778272889dSSatish Balay 
10788272889dSSatish Balay   PetscFunctionBegin;
10798272889dSSatish Balay   if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element");
1080*e6a796c3SToby Isaac   newton = (PetscBool) (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA);
1081*e6a796c3SToby Isaac   ierr = PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, 0., 0., x, w, newton);CHKERRQ(ierr);
10828272889dSSatish Balay   PetscFunctionReturn(0);
10838272889dSSatish Balay }
10848272889dSSatish Balay 
1085744bafbcSMatthew G. Knepley /*@
1086744bafbcSMatthew G. Knepley   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
1087744bafbcSMatthew G. Knepley 
1088744bafbcSMatthew G. Knepley   Not Collective
1089744bafbcSMatthew G. Knepley 
1090744bafbcSMatthew G. Knepley   Input Arguments:
1091744bafbcSMatthew G. Knepley + dim     - The spatial dimension
1092a6b92713SMatthew G. Knepley . Nc      - The number of components
1093744bafbcSMatthew G. Knepley . npoints - number of points in one dimension
1094744bafbcSMatthew G. Knepley . a       - left end of interval (often-1)
1095744bafbcSMatthew G. Knepley - b       - right end of interval (often +1)
1096744bafbcSMatthew G. Knepley 
1097744bafbcSMatthew G. Knepley   Output Argument:
1098744bafbcSMatthew G. Knepley . q - A PetscQuadrature object
1099744bafbcSMatthew G. Knepley 
1100744bafbcSMatthew G. Knepley   Level: intermediate
1101744bafbcSMatthew G. Knepley 
1102744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
1103744bafbcSMatthew G. Knepley @*/
1104a6b92713SMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1105744bafbcSMatthew G. Knepley {
1106a6b92713SMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
1107744bafbcSMatthew G. Knepley   PetscReal     *x, *w, *xw, *ww;
1108744bafbcSMatthew G. Knepley   PetscErrorCode ierr;
1109744bafbcSMatthew G. Knepley 
1110744bafbcSMatthew G. Knepley   PetscFunctionBegin;
1111744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
1112a6b92713SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr);
1113744bafbcSMatthew G. Knepley   /* Set up the Golub-Welsch system */
1114744bafbcSMatthew G. Knepley   switch (dim) {
1115744bafbcSMatthew G. Knepley   case 0:
1116744bafbcSMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
1117744bafbcSMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
1118744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
1119a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
1120744bafbcSMatthew G. Knepley     x[0] = 0.0;
1121a6b92713SMatthew G. Knepley     for (c = 0; c < Nc; ++c) w[c] = 1.0;
1122744bafbcSMatthew G. Knepley     break;
1123744bafbcSMatthew G. Knepley   case 1:
1124a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr);
1125a6b92713SMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr);
1126a6b92713SMatthew G. Knepley     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i];
1127a6b92713SMatthew G. Knepley     ierr = PetscFree(ww);CHKERRQ(ierr);
1128744bafbcSMatthew G. Knepley     break;
1129744bafbcSMatthew G. Knepley   case 2:
1130744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
1131744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
1132744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
1133744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
1134744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+0] = xw[i];
1135744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+1] = xw[j];
1136a6b92713SMatthew G. Knepley         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j];
1137744bafbcSMatthew G. Knepley       }
1138744bafbcSMatthew G. Knepley     }
1139744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
1140744bafbcSMatthew G. Knepley     break;
1141744bafbcSMatthew G. Knepley   case 3:
1142744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
1143744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
1144744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
1145744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
1146744bafbcSMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
1147744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
1148744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
1149744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
1150a6b92713SMatthew G. Knepley           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k];
1151744bafbcSMatthew G. Knepley         }
1152744bafbcSMatthew G. Knepley       }
1153744bafbcSMatthew G. Knepley     }
1154744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
1155744bafbcSMatthew G. Knepley     break;
1156744bafbcSMatthew G. Knepley   default:
1157744bafbcSMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
1158744bafbcSMatthew G. Knepley   }
1159744bafbcSMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
11602f5fb066SToby Isaac   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
1161a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
1162d9bac1caSLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr);
1163744bafbcSMatthew G. Knepley   PetscFunctionReturn(0);
1164744bafbcSMatthew G. Knepley }
1165744bafbcSMatthew G. Knepley 
1166494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
1167494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
1168494e7359SMatthew G. Knepley {
1169494e7359SMatthew G. Knepley   PetscFunctionBegin;
1170494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
1171494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
1172494e7359SMatthew G. Knepley   *zeta = z;
1173494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
1174494e7359SMatthew G. Knepley }
1175494e7359SMatthew G. Knepley 
1176494e7359SMatthew G. Knepley 
1177f5f57ec0SBarry Smith /*@
1178*e6a796c3SToby Isaac   PetscDTStroudConicalQuadrature - create Stroud conical quadrature for a simplex
1179494e7359SMatthew G. Knepley 
1180494e7359SMatthew G. Knepley   Not Collective
1181494e7359SMatthew G. Knepley 
1182494e7359SMatthew G. Knepley   Input Arguments:
1183494e7359SMatthew G. Knepley + dim     - The simplex dimension
1184a6b92713SMatthew G. Knepley . Nc      - The number of components
1185dcce0ee2SMatthew G. Knepley . npoints - The number of points in one dimension
1186494e7359SMatthew G. Knepley . a       - left end of interval (often-1)
1187494e7359SMatthew G. Knepley - b       - right end of interval (often +1)
1188494e7359SMatthew G. Knepley 
1189744bafbcSMatthew G. Knepley   Output Argument:
1190552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
1191494e7359SMatthew G. Knepley 
1192494e7359SMatthew G. Knepley   Level: intermediate
1193494e7359SMatthew G. Knepley 
1194494e7359SMatthew G. Knepley   References:
119596a0c994SBarry Smith .  1. - Karniadakis and Sherwin.  FIAT
1196494e7359SMatthew G. Knepley 
1197*e6a796c3SToby Isaac   Note: For dim == 1, this is Gauss-Legendre quadrature
1198*e6a796c3SToby Isaac 
1199744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
1200494e7359SMatthew G. Knepley @*/
1201*e6a796c3SToby Isaac PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1202494e7359SMatthew G. Knepley {
1203dcce0ee2SMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints;
1204494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
1205*e6a796c3SToby Isaac   PetscInt       i, j, k, c; PetscErrorCode ierr;
1206494e7359SMatthew G. Knepley 
1207494e7359SMatthew G. Knepley   PetscFunctionBegin;
1208494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
1209dcce0ee2SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr);
1210dcce0ee2SMatthew G. Knepley   ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr);
1211494e7359SMatthew G. Knepley   switch (dim) {
1212707aa5c5SMatthew G. Knepley   case 0:
1213707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
1214707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
1215785e854fSJed Brown     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
1216a6b92713SMatthew G. Knepley     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
1217707aa5c5SMatthew G. Knepley     x[0] = 0.0;
1218a6b92713SMatthew G. Knepley     for (c = 0; c < Nc; ++c) w[c] = 1.0;
1219707aa5c5SMatthew G. Knepley     break;
1220494e7359SMatthew G. Knepley   case 1:
1221dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr);
1222*e6a796c3SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr);
1223dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i];
1224a6b92713SMatthew G. Knepley     ierr = PetscFree(wx);CHKERRQ(ierr);
1225494e7359SMatthew G. Knepley     break;
1226494e7359SMatthew G. Knepley   case 2:
1227dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr);
1228*e6a796c3SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
1229*e6a796c3SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
1230dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
1231dcce0ee2SMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
1232dcce0ee2SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr);
1233dcce0ee2SMatthew G. Knepley         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j];
1234494e7359SMatthew G. Knepley       }
1235494e7359SMatthew G. Knepley     }
1236494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
1237494e7359SMatthew G. Knepley     break;
1238494e7359SMatthew G. Knepley   case 3:
1239dcce0ee2SMatthew G. Knepley     ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr);
1240*e6a796c3SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
1241*e6a796c3SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
1242*e6a796c3SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
1243dcce0ee2SMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
1244dcce0ee2SMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
1245dcce0ee2SMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
1246dcce0ee2SMatthew 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);
1247dcce0ee2SMatthew G. Knepley           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k];
1248494e7359SMatthew G. Knepley         }
1249494e7359SMatthew G. Knepley       }
1250494e7359SMatthew G. Knepley     }
1251494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
1252494e7359SMatthew G. Knepley     break;
1253494e7359SMatthew G. Knepley   default:
1254494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
1255494e7359SMatthew G. Knepley   }
125621454ff5SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
12572f5fb066SToby Isaac   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
1258dcce0ee2SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
1259d9bac1caSLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr);
1260494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
1261494e7359SMatthew G. Knepley }
1262494e7359SMatthew G. Knepley 
1263f5f57ec0SBarry Smith /*@
1264b3c0f97bSTom Klotz   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
1265b3c0f97bSTom Klotz 
1266b3c0f97bSTom Klotz   Not Collective
1267b3c0f97bSTom Klotz 
1268b3c0f97bSTom Klotz   Input Arguments:
1269b3c0f97bSTom Klotz + dim   - The cell dimension
1270b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l
1271b3c0f97bSTom Klotz . a     - left end of interval (often-1)
1272b3c0f97bSTom Klotz - b     - right end of interval (often +1)
1273b3c0f97bSTom Klotz 
1274b3c0f97bSTom Klotz   Output Argument:
1275b3c0f97bSTom Klotz . q - A PetscQuadrature object
1276b3c0f97bSTom Klotz 
1277b3c0f97bSTom Klotz   Level: intermediate
1278b3c0f97bSTom Klotz 
1279b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature()
1280b3c0f97bSTom Klotz @*/
1281b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
1282b3c0f97bSTom Klotz {
1283b3c0f97bSTom Klotz   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
1284b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;                  /* Half-width of the integration interval */
1285b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;                  /* Center of the integration interval */
1286b3c0f97bSTom Klotz   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
1287d84b4d08SMatthew G. Knepley   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
1288b3c0f97bSTom Klotz   PetscReal       wk    = 0.5*PETSC_PI;              /* Quadrature weight at x_k */
1289b3c0f97bSTom Klotz   PetscReal      *x, *w;
1290b3c0f97bSTom Klotz   PetscInt        K, k, npoints;
1291b3c0f97bSTom Klotz   PetscErrorCode  ierr;
1292b3c0f97bSTom Klotz 
1293b3c0f97bSTom Klotz   PetscFunctionBegin;
1294b3c0f97bSTom Klotz   if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
1295b3c0f97bSTom Klotz   if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
1296b3c0f97bSTom Klotz   /* Find K such that the weights are < 32 digits of precision */
1297b3c0f97bSTom Klotz   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
12989add2064SThomas Klotz     wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
1299b3c0f97bSTom Klotz   }
1300b3c0f97bSTom Klotz   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
1301b3c0f97bSTom Klotz   ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr);
1302b3c0f97bSTom Klotz   npoints = 2*K-1;
1303b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
1304b3c0f97bSTom Klotz   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
1305b3c0f97bSTom Klotz   /* Center term */
1306b3c0f97bSTom Klotz   x[0] = beta;
1307b3c0f97bSTom Klotz   w[0] = 0.5*alpha*PETSC_PI;
1308b3c0f97bSTom Klotz   for (k = 1; k < K; ++k) {
13099add2064SThomas Klotz     wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
13101118d4bcSLisandro Dalcin     xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h));
1311b3c0f97bSTom Klotz     x[2*k-1] = -alpha*xk+beta;
1312b3c0f97bSTom Klotz     w[2*k-1] = wk;
1313b3c0f97bSTom Klotz     x[2*k+0] =  alpha*xk+beta;
1314b3c0f97bSTom Klotz     w[2*k+0] = wk;
1315b3c0f97bSTom Klotz   }
1316a6b92713SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr);
1317b3c0f97bSTom Klotz   PetscFunctionReturn(0);
1318b3c0f97bSTom Klotz }
1319b3c0f97bSTom Klotz 
1320b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1321b3c0f97bSTom Klotz {
1322b3c0f97bSTom Klotz   const PetscInt  p     = 16;        /* Digits of precision in the evaluation */
1323b3c0f97bSTom Klotz   const PetscReal alpha = (b-a)/2.;  /* Half-width of the integration interval */
1324b3c0f97bSTom Klotz   const PetscReal beta  = (b+a)/2.;  /* Center of the integration interval */
1325b3c0f97bSTom Klotz   PetscReal       h     = 1.0;       /* Step size, length between x_k */
1326b3c0f97bSTom Klotz   PetscInt        l     = 0;         /* Level of refinement, h = 2^{-l} */
1327b3c0f97bSTom Klotz   PetscReal       osum  = 0.0;       /* Integral on last level */
1328b3c0f97bSTom Klotz   PetscReal       psum  = 0.0;       /* Integral on the level before the last level */
1329b3c0f97bSTom Klotz   PetscReal       sum;               /* Integral on current level */
1330446c295cSMatthew G. Knepley   PetscReal       yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1331b3c0f97bSTom Klotz   PetscReal       lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1332b3c0f97bSTom Klotz   PetscReal       wk;                /* Quadrature weight at x_k */
1333b3c0f97bSTom Klotz   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
1334b3c0f97bSTom Klotz   PetscInt        d;                 /* Digits of precision in the integral */
1335b3c0f97bSTom Klotz 
1336b3c0f97bSTom Klotz   PetscFunctionBegin;
1337b3c0f97bSTom Klotz   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1338b3c0f97bSTom Klotz   /* Center term */
1339b3c0f97bSTom Klotz   func(beta, &lval);
1340b3c0f97bSTom Klotz   sum = 0.5*alpha*PETSC_PI*lval;
1341b3c0f97bSTom Klotz   /* */
1342b3c0f97bSTom Klotz   do {
1343b3c0f97bSTom Klotz     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
1344b3c0f97bSTom Klotz     PetscInt  k = 1;
1345b3c0f97bSTom Klotz 
1346b3c0f97bSTom Klotz     ++l;
1347b3c0f97bSTom Klotz     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1348b3c0f97bSTom Klotz     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1349b3c0f97bSTom Klotz     psum = osum;
1350b3c0f97bSTom Klotz     osum = sum;
1351b3c0f97bSTom Klotz     h   *= 0.5;
1352b3c0f97bSTom Klotz     sum *= 0.5;
1353b3c0f97bSTom Klotz     do {
13549add2064SThomas Klotz       wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1355446c295cSMatthew G. Knepley       yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1356446c295cSMatthew G. Knepley       lx = -alpha*(1.0 - yk)+beta;
1357446c295cSMatthew G. Knepley       rx =  alpha*(1.0 - yk)+beta;
1358b3c0f97bSTom Klotz       func(lx, &lval);
1359b3c0f97bSTom Klotz       func(rx, &rval);
1360b3c0f97bSTom Klotz       lterm   = alpha*wk*lval;
1361b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
1362b3c0f97bSTom Klotz       sum    += lterm;
1363b3c0f97bSTom Klotz       rterm   = alpha*wk*rval;
1364b3c0f97bSTom Klotz       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
1365b3c0f97bSTom Klotz       sum    += rterm;
1366b3c0f97bSTom Klotz       ++k;
1367b3c0f97bSTom Klotz       /* Only need to evaluate every other point on refined levels */
1368b3c0f97bSTom Klotz       if (l != 1) ++k;
13699add2064SThomas Klotz     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
1370b3c0f97bSTom Klotz 
1371b3c0f97bSTom Klotz     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
1372b3c0f97bSTom Klotz     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
1373b3c0f97bSTom Klotz     d3 = PetscLog10Real(maxTerm) - p;
137409d48545SBarry Smith     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
137509d48545SBarry Smith     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
1376b3c0f97bSTom Klotz     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
13779add2064SThomas Klotz   } while (d < digits && l < 12);
1378b3c0f97bSTom Klotz   *sol = sum;
1379e510cb1fSThomas Klotz 
1380b3c0f97bSTom Klotz   PetscFunctionReturn(0);
1381b3c0f97bSTom Klotz }
1382b3c0f97bSTom Klotz 
1383497880caSRichard Tran Mills #if defined(PETSC_HAVE_MPFR)
138429f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
138529f144ccSMatthew G. Knepley {
1386e510cb1fSThomas Klotz   const PetscInt  safetyFactor = 2;  /* Calculate abcissa until 2*p digits */
138729f144ccSMatthew G. Knepley   PetscInt        l            = 0;  /* Level of refinement, h = 2^{-l} */
138829f144ccSMatthew G. Knepley   mpfr_t          alpha;             /* Half-width of the integration interval */
138929f144ccSMatthew G. Knepley   mpfr_t          beta;              /* Center of the integration interval */
139029f144ccSMatthew G. Knepley   mpfr_t          h;                 /* Step size, length between x_k */
139129f144ccSMatthew G. Knepley   mpfr_t          osum;              /* Integral on last level */
139229f144ccSMatthew G. Knepley   mpfr_t          psum;              /* Integral on the level before the last level */
139329f144ccSMatthew G. Knepley   mpfr_t          sum;               /* Integral on current level */
139429f144ccSMatthew G. Knepley   mpfr_t          yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
139529f144ccSMatthew G. Knepley   mpfr_t          lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
139629f144ccSMatthew G. Knepley   mpfr_t          wk;                /* Quadrature weight at x_k */
139729f144ccSMatthew G. Knepley   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
139829f144ccSMatthew G. Knepley   PetscInt        d;                 /* Digits of precision in the integral */
139929f144ccSMatthew G. Knepley   mpfr_t          pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
140029f144ccSMatthew G. Knepley 
140129f144ccSMatthew G. Knepley   PetscFunctionBegin;
140229f144ccSMatthew G. Knepley   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
140329f144ccSMatthew G. Knepley   /* Create high precision storage */
1404c9f744b5SMatthew 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);
140529f144ccSMatthew G. Knepley   /* Initialization */
140629f144ccSMatthew G. Knepley   mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
140729f144ccSMatthew G. Knepley   mpfr_set_d(beta,  0.5*(b+a), MPFR_RNDN);
140829f144ccSMatthew G. Knepley   mpfr_set_d(osum,  0.0,       MPFR_RNDN);
140929f144ccSMatthew G. Knepley   mpfr_set_d(psum,  0.0,       MPFR_RNDN);
141029f144ccSMatthew G. Knepley   mpfr_set_d(h,     1.0,       MPFR_RNDN);
141129f144ccSMatthew G. Knepley   mpfr_const_pi(pi2, MPFR_RNDN);
141229f144ccSMatthew G. Knepley   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
141329f144ccSMatthew G. Knepley   /* Center term */
141429f144ccSMatthew G. Knepley   func(0.5*(b+a), &lval);
141529f144ccSMatthew G. Knepley   mpfr_set(sum, pi2, MPFR_RNDN);
141629f144ccSMatthew G. Knepley   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
141729f144ccSMatthew G. Knepley   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
141829f144ccSMatthew G. Knepley   /* */
141929f144ccSMatthew G. Knepley   do {
142029f144ccSMatthew G. Knepley     PetscReal d1, d2, d3, d4;
142129f144ccSMatthew G. Knepley     PetscInt  k = 1;
142229f144ccSMatthew G. Knepley 
142329f144ccSMatthew G. Knepley     ++l;
142429f144ccSMatthew G. Knepley     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
142529f144ccSMatthew G. Knepley     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
142629f144ccSMatthew G. Knepley     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
142729f144ccSMatthew G. Knepley     mpfr_set(psum, osum, MPFR_RNDN);
142829f144ccSMatthew G. Knepley     mpfr_set(osum,  sum, MPFR_RNDN);
142929f144ccSMatthew G. Knepley     mpfr_mul_d(h,   h,   0.5, MPFR_RNDN);
143029f144ccSMatthew G. Knepley     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
143129f144ccSMatthew G. Knepley     do {
143229f144ccSMatthew G. Knepley       mpfr_set_si(kh, k, MPFR_RNDN);
143329f144ccSMatthew G. Knepley       mpfr_mul(kh, kh, h, MPFR_RNDN);
143429f144ccSMatthew G. Knepley       /* Weight */
143529f144ccSMatthew G. Knepley       mpfr_set(wk, h, MPFR_RNDN);
143629f144ccSMatthew G. Knepley       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
143729f144ccSMatthew G. Knepley       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
143829f144ccSMatthew G. Knepley       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
143929f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
144029f144ccSMatthew G. Knepley       mpfr_sqr(tmp, tmp, MPFR_RNDN);
144129f144ccSMatthew G. Knepley       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
144229f144ccSMatthew G. Knepley       mpfr_div(wk, wk, tmp, MPFR_RNDN);
144329f144ccSMatthew G. Knepley       /* Abscissa */
144429f144ccSMatthew G. Knepley       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
144529f144ccSMatthew G. Knepley       mpfr_cosh(tmp, msinh, MPFR_RNDN);
144629f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
144729f144ccSMatthew G. Knepley       mpfr_exp(tmp, msinh, MPFR_RNDN);
144829f144ccSMatthew G. Knepley       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
144929f144ccSMatthew G. Knepley       /* Quadrature points */
145029f144ccSMatthew G. Knepley       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
145129f144ccSMatthew G. Knepley       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
145229f144ccSMatthew G. Knepley       mpfr_add(lx, lx, beta, MPFR_RNDU);
145329f144ccSMatthew G. Knepley       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
145429f144ccSMatthew G. Knepley       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
145529f144ccSMatthew G. Knepley       mpfr_add(rx, rx, beta, MPFR_RNDD);
145629f144ccSMatthew G. Knepley       /* Evaluation */
145729f144ccSMatthew G. Knepley       func(mpfr_get_d(lx, MPFR_RNDU), &lval);
145829f144ccSMatthew G. Knepley       func(mpfr_get_d(rx, MPFR_RNDD), &rval);
145929f144ccSMatthew G. Knepley       /* Update */
146029f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
146129f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
146229f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
146329f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
146429f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
146529f144ccSMatthew G. Knepley       mpfr_set(curTerm, tmp, MPFR_RNDN);
146629f144ccSMatthew G. Knepley       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
146729f144ccSMatthew G. Knepley       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
146829f144ccSMatthew G. Knepley       mpfr_add(sum, sum, tmp, MPFR_RNDN);
146929f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
147029f144ccSMatthew G. Knepley       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
147129f144ccSMatthew G. Knepley       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
147229f144ccSMatthew G. Knepley       ++k;
147329f144ccSMatthew G. Knepley       /* Only need to evaluate every other point on refined levels */
147429f144ccSMatthew G. Knepley       if (l != 1) ++k;
147529f144ccSMatthew G. Knepley       mpfr_log10(tmp, wk, MPFR_RNDN);
147629f144ccSMatthew G. Knepley       mpfr_abs(tmp, tmp, MPFR_RNDN);
1477c9f744b5SMatthew G. Knepley     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
147829f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
147929f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
148029f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
148129f144ccSMatthew G. Knepley     d1 = mpfr_get_d(tmp, MPFR_RNDN);
148229f144ccSMatthew G. Knepley     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
148329f144ccSMatthew G. Knepley     mpfr_abs(tmp, tmp, MPFR_RNDN);
148429f144ccSMatthew G. Knepley     mpfr_log10(tmp, tmp, MPFR_RNDN);
148529f144ccSMatthew G. Knepley     d2 = mpfr_get_d(tmp, MPFR_RNDN);
148629f144ccSMatthew G. Knepley     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
1487c9f744b5SMatthew G. Knepley     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
148829f144ccSMatthew G. Knepley     mpfr_log10(tmp, curTerm, MPFR_RNDN);
148929f144ccSMatthew G. Knepley     d4 = mpfr_get_d(tmp, MPFR_RNDN);
149029f144ccSMatthew G. Knepley     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1491b0649871SThomas Klotz   } while (d < digits && l < 8);
149229f144ccSMatthew G. Knepley   *sol = mpfr_get_d(sum, MPFR_RNDN);
149329f144ccSMatthew G. Knepley   /* Cleanup */
149429f144ccSMatthew G. Knepley   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
149529f144ccSMatthew G. Knepley   PetscFunctionReturn(0);
149629f144ccSMatthew G. Knepley }
1497d525116cSMatthew G. Knepley #else
1498fbfcfee5SBarry Smith 
1499d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1500d525116cSMatthew G. Knepley {
1501d525116cSMatthew G. Knepley   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
1502d525116cSMatthew G. Knepley }
150329f144ccSMatthew G. Knepley #endif
150429f144ccSMatthew G. Knepley 
1505194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
1506194825f6SJed Brown  * A in column-major format
1507194825f6SJed Brown  * Ainv in row-major format
1508194825f6SJed Brown  * tau has length m
1509194825f6SJed Brown  * worksize must be >= max(1,n)
1510194825f6SJed Brown  */
1511194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1512194825f6SJed Brown {
1513194825f6SJed Brown   PetscErrorCode ierr;
1514194825f6SJed Brown   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
1515194825f6SJed Brown   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
1516194825f6SJed Brown 
1517194825f6SJed Brown   PetscFunctionBegin;
1518194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
1519194825f6SJed Brown   {
1520194825f6SJed Brown     PetscInt i,j;
1521dcca6d9dSJed Brown     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
1522194825f6SJed Brown     for (j=0; j<n; j++) {
1523194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1524194825f6SJed Brown     }
1525194825f6SJed Brown     mstride = m;
1526194825f6SJed Brown   }
1527194825f6SJed Brown #else
1528194825f6SJed Brown   A = A_in;
1529194825f6SJed Brown   Ainv = Ainv_out;
1530194825f6SJed Brown #endif
1531194825f6SJed Brown 
1532194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
1533194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
1534194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
1535194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
1536194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1537001a771dSBarry Smith   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1538194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1539194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1540194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1541194825f6SJed Brown 
1542194825f6SJed Brown   /* Extract an explicit representation of Q */
1543194825f6SJed Brown   Q = Ainv;
1544580bdb30SBarry Smith   ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr);
1545194825f6SJed Brown   K = N;                        /* full rank */
1546c964aadfSJose E. Roman   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1547194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1548194825f6SJed Brown 
1549194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1550194825f6SJed Brown   Alpha = 1.0;
1551194825f6SJed Brown   ldb = lda;
1552001a771dSBarry Smith   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1553194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
1554194825f6SJed Brown 
1555194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
1556194825f6SJed Brown   {
1557194825f6SJed Brown     PetscInt i;
1558194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1559194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
1560194825f6SJed Brown   }
1561194825f6SJed Brown #endif
1562194825f6SJed Brown   PetscFunctionReturn(0);
1563194825f6SJed Brown }
1564194825f6SJed Brown 
1565194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1566194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1567194825f6SJed Brown {
1568194825f6SJed Brown   PetscErrorCode ierr;
1569194825f6SJed Brown   PetscReal      *Bv;
1570194825f6SJed Brown   PetscInt       i,j;
1571194825f6SJed Brown 
1572194825f6SJed Brown   PetscFunctionBegin;
1573785e854fSJed Brown   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
1574194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
1575194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
1576194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1577194825f6SJed Brown   for (i=0; i<ninterval; i++) {
1578194825f6SJed Brown     for (j=0; j<ndegree; j++) {
1579194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1580194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1581194825f6SJed Brown     }
1582194825f6SJed Brown   }
1583194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
1584194825f6SJed Brown   PetscFunctionReturn(0);
1585194825f6SJed Brown }
1586194825f6SJed Brown 
1587194825f6SJed Brown /*@
1588194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
1589194825f6SJed Brown 
1590194825f6SJed Brown    Not Collective
1591194825f6SJed Brown 
1592194825f6SJed Brown    Input Arguments:
1593194825f6SJed Brown +  degree - degree of reconstruction polynomial
1594194825f6SJed Brown .  nsource - number of source intervals
1595194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1596194825f6SJed Brown .  ntarget - number of target intervals
1597194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
1598194825f6SJed Brown 
1599194825f6SJed Brown    Output Arguments:
1600194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
1601194825f6SJed Brown 
1602194825f6SJed Brown    Level: advanced
1603194825f6SJed Brown 
1604194825f6SJed Brown .seealso: PetscDTLegendreEval()
1605194825f6SJed Brown @*/
1606194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1607194825f6SJed Brown {
1608194825f6SJed Brown   PetscErrorCode ierr;
1609194825f6SJed Brown   PetscInt       i,j,k,*bdegrees,worksize;
1610194825f6SJed Brown   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1611194825f6SJed Brown   PetscScalar    *tau,*work;
1612194825f6SJed Brown 
1613194825f6SJed Brown   PetscFunctionBegin;
1614194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
1615194825f6SJed Brown   PetscValidRealPointer(targetx,5);
1616194825f6SJed Brown   PetscValidRealPointer(R,6);
1617194825f6SJed 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);
1618194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
1619194825f6SJed Brown   for (i=0; i<nsource; i++) {
162057622a8eSBarry 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]);
1621194825f6SJed Brown   }
1622194825f6SJed Brown   for (i=0; i<ntarget; i++) {
162357622a8eSBarry 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]);
1624194825f6SJed Brown   }
1625194825f6SJed Brown #endif
1626194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
1627194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1628194825f6SJed Brown   center = (xmin + xmax)/2;
1629194825f6SJed Brown   hscale = (xmax - xmin)/2;
1630194825f6SJed Brown   worksize = nsource;
1631dcca6d9dSJed Brown   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
1632dcca6d9dSJed Brown   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
1633194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1634194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1635194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
1636194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
1637194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1638194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
1639194825f6SJed Brown   for (i=0; i<ntarget; i++) {
1640194825f6SJed Brown     PetscReal rowsum = 0;
1641194825f6SJed Brown     for (j=0; j<nsource; j++) {
1642194825f6SJed Brown       PetscReal sum = 0;
1643194825f6SJed Brown       for (k=0; k<degree+1; k++) {
1644194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1645194825f6SJed Brown       }
1646194825f6SJed Brown       R[i*nsource+j] = sum;
1647194825f6SJed Brown       rowsum += sum;
1648194825f6SJed Brown     }
1649194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1650194825f6SJed Brown   }
1651194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
1652194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
1653194825f6SJed Brown   PetscFunctionReturn(0);
1654194825f6SJed Brown }
1655916e780bShannah_mairs 
1656916e780bShannah_mairs /*@C
1657916e780bShannah_mairs    PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points
1658916e780bShannah_mairs 
1659916e780bShannah_mairs    Not Collective
1660916e780bShannah_mairs 
1661916e780bShannah_mairs    Input Parameter:
1662916e780bShannah_mairs +  n - the number of GLL nodes
1663916e780bShannah_mairs .  nodes - the GLL nodes
1664916e780bShannah_mairs .  weights - the GLL weights
1665916e780bShannah_mairs .  f - the function values at the nodes
1666916e780bShannah_mairs 
1667916e780bShannah_mairs    Output Parameter:
1668916e780bShannah_mairs .  in - the value of the integral
1669916e780bShannah_mairs 
1670916e780bShannah_mairs    Level: beginner
1671916e780bShannah_mairs 
1672916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature()
1673916e780bShannah_mairs 
1674916e780bShannah_mairs @*/
1675916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in)
1676916e780bShannah_mairs {
1677916e780bShannah_mairs   PetscInt          i;
1678916e780bShannah_mairs 
1679916e780bShannah_mairs   PetscFunctionBegin;
1680916e780bShannah_mairs   *in = 0.;
1681916e780bShannah_mairs   for (i=0; i<n; i++) {
1682916e780bShannah_mairs     *in += f[i]*f[i]*weights[i];
1683916e780bShannah_mairs   }
1684916e780bShannah_mairs   PetscFunctionReturn(0);
1685916e780bShannah_mairs }
1686916e780bShannah_mairs 
1687916e780bShannah_mairs /*@C
1688916e780bShannah_mairs    PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element
1689916e780bShannah_mairs 
1690916e780bShannah_mairs    Not Collective
1691916e780bShannah_mairs 
1692916e780bShannah_mairs    Input Parameter:
1693916e780bShannah_mairs +  n - the number of GLL nodes
1694916e780bShannah_mairs .  nodes - the GLL nodes
1695916e780bShannah_mairs .  weights - the GLL weights
1696916e780bShannah_mairs 
1697916e780bShannah_mairs    Output Parameter:
1698916e780bShannah_mairs .  A - the stiffness element
1699916e780bShannah_mairs 
1700916e780bShannah_mairs    Level: beginner
1701916e780bShannah_mairs 
1702916e780bShannah_mairs    Notes:
1703916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy()
1704916e780bShannah_mairs 
1705916e780bShannah_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)
1706916e780bShannah_mairs 
1707916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1708916e780bShannah_mairs 
1709916e780bShannah_mairs @*/
1710916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1711916e780bShannah_mairs {
1712916e780bShannah_mairs   PetscReal        **A;
1713916e780bShannah_mairs   PetscErrorCode  ierr;
1714916e780bShannah_mairs   const PetscReal  *gllnodes = nodes;
1715916e780bShannah_mairs   const PetscInt   p = n-1;
1716916e780bShannah_mairs   PetscReal        z0,z1,z2 = -1,x,Lpj,Lpr;
1717916e780bShannah_mairs   PetscInt         i,j,nn,r;
1718916e780bShannah_mairs 
1719916e780bShannah_mairs   PetscFunctionBegin;
1720916e780bShannah_mairs   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1721916e780bShannah_mairs   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1722916e780bShannah_mairs   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1723916e780bShannah_mairs 
1724916e780bShannah_mairs   for (j=1; j<p; j++) {
1725916e780bShannah_mairs     x  = gllnodes[j];
1726916e780bShannah_mairs     z0 = 1.;
1727916e780bShannah_mairs     z1 = x;
1728916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1729916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1730916e780bShannah_mairs       z0 = z1;
1731916e780bShannah_mairs       z1 = z2;
1732916e780bShannah_mairs     }
1733916e780bShannah_mairs     Lpj=z2;
1734916e780bShannah_mairs     for (r=1; r<p; r++) {
1735916e780bShannah_mairs       if (r == j) {
1736916e780bShannah_mairs         A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj);
1737916e780bShannah_mairs       } else {
1738916e780bShannah_mairs         x  = gllnodes[r];
1739916e780bShannah_mairs         z0 = 1.;
1740916e780bShannah_mairs         z1 = x;
1741916e780bShannah_mairs         for (nn=1; nn<p; nn++) {
1742916e780bShannah_mairs           z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1743916e780bShannah_mairs           z0 = z1;
1744916e780bShannah_mairs           z1 = z2;
1745916e780bShannah_mairs         }
1746916e780bShannah_mairs         Lpr     = z2;
1747916e780bShannah_mairs         A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r]));
1748916e780bShannah_mairs       }
1749916e780bShannah_mairs     }
1750916e780bShannah_mairs   }
1751916e780bShannah_mairs   for (j=1; j<p+1; j++) {
1752916e780bShannah_mairs     x  = gllnodes[j];
1753916e780bShannah_mairs     z0 = 1.;
1754916e780bShannah_mairs     z1 = x;
1755916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1756916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1757916e780bShannah_mairs       z0 = z1;
1758916e780bShannah_mairs       z1 = z2;
1759916e780bShannah_mairs     }
1760916e780bShannah_mairs     Lpj     = z2;
1761916e780bShannah_mairs     A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j]));
1762916e780bShannah_mairs     A[0][j] = A[j][0];
1763916e780bShannah_mairs   }
1764916e780bShannah_mairs   for (j=0; j<p; j++) {
1765916e780bShannah_mairs     x  = gllnodes[j];
1766916e780bShannah_mairs     z0 = 1.;
1767916e780bShannah_mairs     z1 = x;
1768916e780bShannah_mairs     for (nn=1; nn<p; nn++) {
1769916e780bShannah_mairs       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1770916e780bShannah_mairs       z0 = z1;
1771916e780bShannah_mairs       z1 = z2;
1772916e780bShannah_mairs     }
1773916e780bShannah_mairs     Lpj=z2;
1774916e780bShannah_mairs 
1775916e780bShannah_mairs     A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j]));
1776916e780bShannah_mairs     A[j][p] = A[p][j];
1777916e780bShannah_mairs   }
1778916e780bShannah_mairs   A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.;
1779916e780bShannah_mairs   A[p][p]=A[0][0];
1780916e780bShannah_mairs   *AA = A;
1781916e780bShannah_mairs   PetscFunctionReturn(0);
1782916e780bShannah_mairs }
1783916e780bShannah_mairs 
1784916e780bShannah_mairs /*@C
1785916e780bShannah_mairs    PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element
1786916e780bShannah_mairs 
1787916e780bShannah_mairs    Not Collective
1788916e780bShannah_mairs 
1789916e780bShannah_mairs    Input Parameter:
1790916e780bShannah_mairs +  n - the number of GLL nodes
1791916e780bShannah_mairs .  nodes - the GLL nodes
1792916e780bShannah_mairs .  weights - the GLL weightss
1793916e780bShannah_mairs -  A - the stiffness element
1794916e780bShannah_mairs 
1795916e780bShannah_mairs    Level: beginner
1796916e780bShannah_mairs 
1797916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate()
1798916e780bShannah_mairs 
1799916e780bShannah_mairs @*/
1800916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1801916e780bShannah_mairs {
1802916e780bShannah_mairs   PetscErrorCode ierr;
1803916e780bShannah_mairs 
1804916e780bShannah_mairs   PetscFunctionBegin;
1805916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1806916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1807916e780bShannah_mairs   *AA  = NULL;
1808916e780bShannah_mairs   PetscFunctionReturn(0);
1809916e780bShannah_mairs }
1810916e780bShannah_mairs 
1811916e780bShannah_mairs /*@C
1812916e780bShannah_mairs    PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element
1813916e780bShannah_mairs 
1814916e780bShannah_mairs    Not Collective
1815916e780bShannah_mairs 
1816916e780bShannah_mairs    Input Parameter:
1817916e780bShannah_mairs +  n - the number of GLL nodes
1818916e780bShannah_mairs .  nodes - the GLL nodes
1819916e780bShannah_mairs .  weights - the GLL weights
1820916e780bShannah_mairs 
1821916e780bShannah_mairs    Output Parameter:
1822916e780bShannah_mairs .  AA - the stiffness element
1823916e780bShannah_mairs -  AAT - the transpose of AA (pass in NULL if you do not need this array)
1824916e780bShannah_mairs 
1825916e780bShannah_mairs    Level: beginner
1826916e780bShannah_mairs 
1827916e780bShannah_mairs    Notes:
1828916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementGradientDestroy()
1829916e780bShannah_mairs 
1830916e780bShannah_mairs    You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1831916e780bShannah_mairs 
1832916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1833916e780bShannah_mairs 
1834916e780bShannah_mairs @*/
1835916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
1836916e780bShannah_mairs {
1837916e780bShannah_mairs   PetscReal        **A, **AT = NULL;
1838916e780bShannah_mairs   PetscErrorCode  ierr;
1839916e780bShannah_mairs   const PetscReal  *gllnodes = nodes;
1840916e780bShannah_mairs   const PetscInt   p = n-1;
1841*e6a796c3SToby Isaac   PetscReal        Li, Lj,d0;
1842916e780bShannah_mairs   PetscInt         i,j;
1843916e780bShannah_mairs 
1844916e780bShannah_mairs   PetscFunctionBegin;
1845916e780bShannah_mairs   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1846916e780bShannah_mairs   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1847916e780bShannah_mairs   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1848916e780bShannah_mairs 
1849916e780bShannah_mairs   if (AAT) {
1850916e780bShannah_mairs     ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr);
1851916e780bShannah_mairs     ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr);
1852916e780bShannah_mairs     for (i=1; i<n; i++) AT[i] = AT[i-1]+n;
1853916e780bShannah_mairs   }
1854916e780bShannah_mairs 
1855916e780bShannah_mairs   if (n==1) {A[0][0] = 0.;}
1856916e780bShannah_mairs   d0 = (PetscReal)p*((PetscReal)p+1.)/4.;
1857916e780bShannah_mairs   for  (i=0; i<n; i++) {
1858916e780bShannah_mairs     for  (j=0; j<n; j++) {
1859916e780bShannah_mairs       A[i][j] = 0.;
1860*e6a796c3SToby Isaac       ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[i], &Li);CHKERRQ(ierr);
1861*e6a796c3SToby Isaac       ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[j], &Lj);CHKERRQ(ierr);
1862916e780bShannah_mairs       if (i!=j)             A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j]));
1863916e780bShannah_mairs       if ((j==i) && (i==0)) A[i][j] = -d0;
1864916e780bShannah_mairs       if (j==i && i==p)     A[i][j] = d0;
1865916e780bShannah_mairs       if (AT) AT[j][i] = A[i][j];
1866916e780bShannah_mairs     }
1867916e780bShannah_mairs   }
1868916e780bShannah_mairs   if (AAT) *AAT = AT;
1869916e780bShannah_mairs   *AA  = A;
1870916e780bShannah_mairs   PetscFunctionReturn(0);
1871916e780bShannah_mairs }
1872916e780bShannah_mairs 
1873916e780bShannah_mairs /*@C
1874916e780bShannah_mairs    PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate()
1875916e780bShannah_mairs 
1876916e780bShannah_mairs    Not Collective
1877916e780bShannah_mairs 
1878916e780bShannah_mairs    Input Parameter:
1879916e780bShannah_mairs +  n - the number of GLL nodes
1880916e780bShannah_mairs .  nodes - the GLL nodes
1881916e780bShannah_mairs .  weights - the GLL weights
1882916e780bShannah_mairs .  AA - the stiffness element
1883916e780bShannah_mairs -  AAT - the transpose of the element
1884916e780bShannah_mairs 
1885916e780bShannah_mairs    Level: beginner
1886916e780bShannah_mairs 
1887916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate()
1888916e780bShannah_mairs 
1889916e780bShannah_mairs @*/
1890916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
1891916e780bShannah_mairs {
1892916e780bShannah_mairs   PetscErrorCode ierr;
1893916e780bShannah_mairs 
1894916e780bShannah_mairs   PetscFunctionBegin;
1895916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1896916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1897916e780bShannah_mairs   *AA  = NULL;
1898916e780bShannah_mairs   if (*AAT) {
1899916e780bShannah_mairs     ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr);
1900916e780bShannah_mairs     ierr = PetscFree(*AAT);CHKERRQ(ierr);
1901916e780bShannah_mairs     *AAT  = NULL;
1902916e780bShannah_mairs   }
1903916e780bShannah_mairs   PetscFunctionReturn(0);
1904916e780bShannah_mairs }
1905916e780bShannah_mairs 
1906916e780bShannah_mairs /*@C
1907916e780bShannah_mairs    PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element
1908916e780bShannah_mairs 
1909916e780bShannah_mairs    Not Collective
1910916e780bShannah_mairs 
1911916e780bShannah_mairs    Input Parameter:
1912916e780bShannah_mairs +  n - the number of GLL nodes
1913916e780bShannah_mairs .  nodes - the GLL nodes
1914916e780bShannah_mairs .  weights - the GLL weightss
1915916e780bShannah_mairs 
1916916e780bShannah_mairs    Output Parameter:
1917916e780bShannah_mairs .  AA - the stiffness element
1918916e780bShannah_mairs 
1919916e780bShannah_mairs    Level: beginner
1920916e780bShannah_mairs 
1921916e780bShannah_mairs    Notes:
1922916e780bShannah_mairs     Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy()
1923916e780bShannah_mairs 
1924916e780bShannah_mairs    This is the same as the Gradient operator multiplied by the diagonal mass matrix
1925916e780bShannah_mairs 
1926916e780bShannah_mairs    You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1927916e780bShannah_mairs 
1928916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy()
1929916e780bShannah_mairs 
1930916e780bShannah_mairs @*/
1931916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1932916e780bShannah_mairs {
1933916e780bShannah_mairs   PetscReal       **D;
1934916e780bShannah_mairs   PetscErrorCode  ierr;
1935916e780bShannah_mairs   const PetscReal  *gllweights = weights;
1936916e780bShannah_mairs   const PetscInt   glln = n;
1937916e780bShannah_mairs   PetscInt         i,j;
1938916e780bShannah_mairs 
1939916e780bShannah_mairs   PetscFunctionBegin;
1940916e780bShannah_mairs   ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr);
1941916e780bShannah_mairs   for (i=0; i<glln; i++){
1942916e780bShannah_mairs     for (j=0; j<glln; j++) {
1943916e780bShannah_mairs       D[i][j] = gllweights[i]*D[i][j];
1944916e780bShannah_mairs     }
1945916e780bShannah_mairs   }
1946916e780bShannah_mairs   *AA = D;
1947916e780bShannah_mairs   PetscFunctionReturn(0);
1948916e780bShannah_mairs }
1949916e780bShannah_mairs 
1950916e780bShannah_mairs /*@C
1951916e780bShannah_mairs    PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element
1952916e780bShannah_mairs 
1953916e780bShannah_mairs    Not Collective
1954916e780bShannah_mairs 
1955916e780bShannah_mairs    Input Parameter:
1956916e780bShannah_mairs +  n - the number of GLL nodes
1957916e780bShannah_mairs .  nodes - the GLL nodes
1958916e780bShannah_mairs .  weights - the GLL weights
1959916e780bShannah_mairs -  A - advection
1960916e780bShannah_mairs 
1961916e780bShannah_mairs    Level: beginner
1962916e780bShannah_mairs 
1963916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate()
1964916e780bShannah_mairs 
1965916e780bShannah_mairs @*/
1966916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1967916e780bShannah_mairs {
1968916e780bShannah_mairs   PetscErrorCode ierr;
1969916e780bShannah_mairs 
1970916e780bShannah_mairs   PetscFunctionBegin;
1971916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1972916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
1973916e780bShannah_mairs   *AA  = NULL;
1974916e780bShannah_mairs   PetscFunctionReturn(0);
1975916e780bShannah_mairs }
1976916e780bShannah_mairs 
1977916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1978916e780bShannah_mairs {
1979916e780bShannah_mairs   PetscReal        **A;
1980916e780bShannah_mairs   PetscErrorCode  ierr;
1981916e780bShannah_mairs   const PetscReal  *gllweights = weights;
1982916e780bShannah_mairs   const PetscInt   glln = n;
1983916e780bShannah_mairs   PetscInt         i,j;
1984916e780bShannah_mairs 
1985916e780bShannah_mairs   PetscFunctionBegin;
1986916e780bShannah_mairs   ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr);
1987916e780bShannah_mairs   ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr);
1988916e780bShannah_mairs   for (i=1; i<glln; i++) A[i] = A[i-1]+glln;
1989916e780bShannah_mairs   if (glln==1) {A[0][0] = 0.;}
1990916e780bShannah_mairs   for  (i=0; i<glln; i++) {
1991916e780bShannah_mairs     for  (j=0; j<glln; j++) {
1992916e780bShannah_mairs       A[i][j] = 0.;
1993916e780bShannah_mairs       if (j==i)     A[i][j] = gllweights[i];
1994916e780bShannah_mairs     }
1995916e780bShannah_mairs   }
1996916e780bShannah_mairs   *AA  = A;
1997916e780bShannah_mairs   PetscFunctionReturn(0);
1998916e780bShannah_mairs }
1999916e780bShannah_mairs 
2000916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2001916e780bShannah_mairs {
2002916e780bShannah_mairs   PetscErrorCode ierr;
2003916e780bShannah_mairs 
2004916e780bShannah_mairs   PetscFunctionBegin;
2005916e780bShannah_mairs   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
2006916e780bShannah_mairs   ierr = PetscFree(*AA);CHKERRQ(ierr);
2007916e780bShannah_mairs   *AA  = NULL;
2008916e780bShannah_mairs   PetscFunctionReturn(0);
2009916e780bShannah_mairs }
2010