xref: /petsc/src/dm/dt/interface/dt.c (revision 40d8ff71d7e79bf1110739faedbf1fcb1334195b)
137045ce4SJed Brown /* Discretization tools */
237045ce4SJed Brown 
3a6fc04d9SSatish Balay #include <petscconf.h>
4a6fc04d9SSatish Balay #if defined(PETSC_HAVE_MATHIMF_H)
5a6fc04d9SSatish Balay #include <mathimf.h>           /* this needs to be included before math.h */
6a6fc04d9SSatish Balay #endif
7a6fc04d9SSatish Balay 
80c35b76eSJed Brown #include <petscdt.h>            /*I "petscdt.h" I*/
937045ce4SJed Brown #include <petscblaslapack.h>
10194825f6SJed Brown #include <petsc-private/petscimpl.h>
1121454ff5SMatthew G. Knepley #include <petsc-private/dtimpl.h>
12665c2dedSJed Brown #include <petscviewer.h>
1359804f93SMatthew G. Knepley #include <petscdmplex.h>
1459804f93SMatthew G. Knepley #include <petscdmshell.h>
1537045ce4SJed Brown 
160bfcf5a5SMatthew G. Knepley static PetscBool GaussCite       = PETSC_FALSE;
170bfcf5a5SMatthew G. Knepley const char       GaussCitation[] = "@article{GolubWelsch1969,\n"
180bfcf5a5SMatthew G. Knepley                                    "  author  = {Golub and Welsch},\n"
190bfcf5a5SMatthew G. Knepley                                    "  title   = {Calculation of Quadrature Rules},\n"
200bfcf5a5SMatthew G. Knepley                                    "  journal = {Math. Comp.},\n"
210bfcf5a5SMatthew G. Knepley                                    "  volume  = {23},\n"
220bfcf5a5SMatthew G. Knepley                                    "  number  = {106},\n"
230bfcf5a5SMatthew G. Knepley                                    "  pages   = {221--230},\n"
240bfcf5a5SMatthew G. Knepley                                    "  year    = {1969}\n}\n";
250bfcf5a5SMatthew G. Knepley 
2637045ce4SJed Brown #undef __FUNCT__
2721454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureCreate"
28*40d8ff71SMatthew G. Knepley /*@
29*40d8ff71SMatthew G. Knepley   PetscQuadratureCreate - Create a PetscQuadrature object
30*40d8ff71SMatthew G. Knepley 
31*40d8ff71SMatthew G. Knepley   Collective on MPI_Comm
32*40d8ff71SMatthew G. Knepley 
33*40d8ff71SMatthew G. Knepley   Input Parameter:
34*40d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object
35*40d8ff71SMatthew G. Knepley 
36*40d8ff71SMatthew G. Knepley   Output Parameter:
37*40d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
38*40d8ff71SMatthew G. Knepley 
39*40d8ff71SMatthew G. Knepley   Level: beginner
40*40d8ff71SMatthew G. Knepley 
41*40d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, create
42*40d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
43*40d8ff71SMatthew G. Knepley @*/
4421454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
4521454ff5SMatthew G. Knepley {
4621454ff5SMatthew G. Knepley   PetscErrorCode ierr;
4721454ff5SMatthew G. Knepley 
4821454ff5SMatthew G. Knepley   PetscFunctionBegin;
4921454ff5SMatthew G. Knepley   PetscValidPointer(q, 2);
5021454ff5SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
5121454ff5SMatthew G. Knepley   ierr = PetscHeaderCreate(*q,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
5221454ff5SMatthew G. Knepley   (*q)->dim       = -1;
5321454ff5SMatthew G. Knepley   (*q)->numPoints = 0;
5421454ff5SMatthew G. Knepley   (*q)->points    = NULL;
5521454ff5SMatthew G. Knepley   (*q)->weights   = NULL;
5621454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
5721454ff5SMatthew G. Knepley }
5821454ff5SMatthew G. Knepley 
5921454ff5SMatthew G. Knepley #undef __FUNCT__
60bfa639d9SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDestroy"
61*40d8ff71SMatthew G. Knepley /*@
62*40d8ff71SMatthew G. Knepley   PetscQuadratureDestroy - Destroys a PetscQuadrature object
63*40d8ff71SMatthew G. Knepley 
64*40d8ff71SMatthew G. Knepley   Collective on PetscQuadrature
65*40d8ff71SMatthew G. Knepley 
66*40d8ff71SMatthew G. Knepley   Input Parameter:
67*40d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
68*40d8ff71SMatthew G. Knepley 
69*40d8ff71SMatthew G. Knepley   Level: beginner
70*40d8ff71SMatthew G. Knepley 
71*40d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, destroy
72*40d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
73*40d8ff71SMatthew G. Knepley @*/
74bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
75bfa639d9SMatthew G. Knepley {
76bfa639d9SMatthew G. Knepley   PetscErrorCode ierr;
77bfa639d9SMatthew G. Knepley 
78bfa639d9SMatthew G. Knepley   PetscFunctionBegin;
7921454ff5SMatthew G. Knepley   if (!*q) PetscFunctionReturn(0);
8021454ff5SMatthew G. Knepley   PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1);
8121454ff5SMatthew G. Knepley   if (--((PetscObject)(*q))->refct > 0) {
8221454ff5SMatthew G. Knepley     *q = NULL;
8321454ff5SMatthew G. Knepley     PetscFunctionReturn(0);
8421454ff5SMatthew G. Knepley   }
8521454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
8621454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
8721454ff5SMatthew G. Knepley   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
8821454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
8921454ff5SMatthew G. Knepley }
9021454ff5SMatthew G. Knepley 
9121454ff5SMatthew G. Knepley #undef __FUNCT__
9221454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetData"
93*40d8ff71SMatthew G. Knepley /*@C
94*40d8ff71SMatthew G. Knepley   PetscQuadratureGetData - Returns the data defining the quadrature
95*40d8ff71SMatthew G. Knepley 
96*40d8ff71SMatthew G. Knepley   Not collective
97*40d8ff71SMatthew G. Knepley 
98*40d8ff71SMatthew G. Knepley   Input Parameter:
99*40d8ff71SMatthew G. Knepley . q  - The PetscQuadrature object
100*40d8ff71SMatthew G. Knepley 
101*40d8ff71SMatthew G. Knepley   Output Parameters:
102*40d8ff71SMatthew G. Knepley + dim - The spatial dimension
103*40d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
104*40d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
105*40d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
106*40d8ff71SMatthew G. Knepley 
107*40d8ff71SMatthew G. Knepley   Level: intermediate
108*40d8ff71SMatthew G. Knepley 
109*40d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature
110*40d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
111*40d8ff71SMatthew G. Knepley @*/
11221454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
11321454ff5SMatthew G. Knepley {
11421454ff5SMatthew G. Knepley   PetscFunctionBegin;
11521454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
11621454ff5SMatthew G. Knepley   if (dim) {
11721454ff5SMatthew G. Knepley     PetscValidPointer(dim, 2);
11821454ff5SMatthew G. Knepley     *dim = q->dim;
11921454ff5SMatthew G. Knepley   }
12021454ff5SMatthew G. Knepley   if (npoints) {
12121454ff5SMatthew G. Knepley     PetscValidPointer(npoints, 3);
12221454ff5SMatthew G. Knepley     *npoints = q->numPoints;
12321454ff5SMatthew G. Knepley   }
12421454ff5SMatthew G. Knepley   if (points) {
12521454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
12621454ff5SMatthew G. Knepley     *points = q->points;
12721454ff5SMatthew G. Knepley   }
12821454ff5SMatthew G. Knepley   if (weights) {
12921454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
13021454ff5SMatthew G. Knepley     *weights = q->weights;
13121454ff5SMatthew G. Knepley   }
13221454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
13321454ff5SMatthew G. Knepley }
13421454ff5SMatthew G. Knepley 
13521454ff5SMatthew G. Knepley #undef __FUNCT__
13621454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetData"
137*40d8ff71SMatthew G. Knepley /*@C
138*40d8ff71SMatthew G. Knepley   PetscQuadratureSetData - Sets the data defining the quadrature
139*40d8ff71SMatthew G. Knepley 
140*40d8ff71SMatthew G. Knepley   Not collective
141*40d8ff71SMatthew G. Knepley 
142*40d8ff71SMatthew G. Knepley   Input Parameters:
143*40d8ff71SMatthew G. Knepley + q  - The PetscQuadrature object
144*40d8ff71SMatthew G. Knepley . dim - The spatial dimension
145*40d8ff71SMatthew G. Knepley . npoints - The number of quadrature points
146*40d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point
147*40d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point
148*40d8ff71SMatthew G. Knepley 
149*40d8ff71SMatthew G. Knepley   Level: intermediate
150*40d8ff71SMatthew G. Knepley 
151*40d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature
152*40d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
153*40d8ff71SMatthew G. Knepley @*/
15421454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
15521454ff5SMatthew G. Knepley {
15621454ff5SMatthew G. Knepley   PetscFunctionBegin;
15721454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
15821454ff5SMatthew G. Knepley   if (dim >= 0)     q->dim       = dim;
15921454ff5SMatthew G. Knepley   if (npoints >= 0) q->numPoints = npoints;
16021454ff5SMatthew G. Knepley   if (points) {
16121454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
16221454ff5SMatthew G. Knepley     q->points = points;
16321454ff5SMatthew G. Knepley   }
16421454ff5SMatthew G. Knepley   if (weights) {
16521454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
16621454ff5SMatthew G. Knepley     q->weights = weights;
16721454ff5SMatthew G. Knepley   }
168f9fd7fdbSMatthew G. Knepley   PetscFunctionReturn(0);
169f9fd7fdbSMatthew G. Knepley }
170f9fd7fdbSMatthew G. Knepley 
171f9fd7fdbSMatthew G. Knepley #undef __FUNCT__
172f9fd7fdbSMatthew G. Knepley #define __FUNCT__ "PetscQuadratureView"
173*40d8ff71SMatthew G. Knepley /*@C
174*40d8ff71SMatthew G. Knepley   PetscQuadratureView - Views a PetscQuadrature object
175*40d8ff71SMatthew G. Knepley 
176*40d8ff71SMatthew G. Knepley   Collective on PetscQuadrature
177*40d8ff71SMatthew G. Knepley 
178*40d8ff71SMatthew G. Knepley   Input Parameters:
179*40d8ff71SMatthew G. Knepley + q  - The PetscQuadrature object
180*40d8ff71SMatthew G. Knepley - viewer - The PetscViewer object
181*40d8ff71SMatthew G. Knepley 
182*40d8ff71SMatthew G. Knepley   Level: beginner
183*40d8ff71SMatthew G. Knepley 
184*40d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, view
185*40d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
186*40d8ff71SMatthew G. Knepley @*/
187f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
188f9fd7fdbSMatthew G. Knepley {
189f9fd7fdbSMatthew G. Knepley   PetscInt       q, d;
190f9fd7fdbSMatthew G. Knepley   PetscErrorCode ierr;
191f9fd7fdbSMatthew G. Knepley 
192f9fd7fdbSMatthew G. Knepley   PetscFunctionBegin;
19398c3331eSBarry Smith   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr);
19421454ff5SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n  (", quad->numPoints);CHKERRQ(ierr);
19521454ff5SMatthew G. Knepley   for (q = 0; q < quad->numPoints; ++q) {
19621454ff5SMatthew G. Knepley     for (d = 0; d < quad->dim; ++d) {
197f9fd7fdbSMatthew G. Knepley       if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);
198ab15ae43SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
199f9fd7fdbSMatthew G. Knepley     }
200ab15ae43SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr);
201f9fd7fdbSMatthew G. Knepley   }
202bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
203bfa639d9SMatthew G. Knepley }
204bfa639d9SMatthew G. Knepley 
205bfa639d9SMatthew G. Knepley #undef __FUNCT__
20637045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval"
20737045ce4SJed Brown /*@
20837045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
20937045ce4SJed Brown 
21037045ce4SJed Brown    Not Collective
21137045ce4SJed Brown 
21237045ce4SJed Brown    Input Arguments:
21337045ce4SJed Brown +  npoints - number of spatial points to evaluate at
21437045ce4SJed Brown .  points - array of locations to evaluate at
21537045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
21637045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
21737045ce4SJed Brown 
21837045ce4SJed Brown    Output Arguments:
2190298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
2200298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
2210298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
22237045ce4SJed Brown 
22337045ce4SJed Brown    Level: intermediate
22437045ce4SJed Brown 
22537045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
22637045ce4SJed Brown @*/
22737045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
22837045ce4SJed Brown {
22937045ce4SJed Brown   PetscInt i,maxdegree;
23037045ce4SJed Brown 
23137045ce4SJed Brown   PetscFunctionBegin;
23237045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
23337045ce4SJed Brown   maxdegree = degrees[ndegree-1];
23437045ce4SJed Brown   for (i=0; i<npoints; i++) {
23537045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
23637045ce4SJed Brown     PetscInt  j,k;
23737045ce4SJed Brown     x    = points[i];
23837045ce4SJed Brown     pm2  = 0;
23937045ce4SJed Brown     pm1  = 1;
24037045ce4SJed Brown     pd2  = 0;
24137045ce4SJed Brown     pd1  = 0;
24237045ce4SJed Brown     pdd2 = 0;
24337045ce4SJed Brown     pdd1 = 0;
24437045ce4SJed Brown     k    = 0;
24537045ce4SJed Brown     if (degrees[k] == 0) {
24637045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
24737045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
24837045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
24937045ce4SJed Brown       k++;
25037045ce4SJed Brown     }
25137045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
25237045ce4SJed Brown       PetscReal p,d,dd;
25337045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
25437045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
25537045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
25637045ce4SJed Brown       pm2  = pm1;
25737045ce4SJed Brown       pm1  = p;
25837045ce4SJed Brown       pd2  = pd1;
25937045ce4SJed Brown       pd1  = d;
26037045ce4SJed Brown       pdd2 = pdd1;
26137045ce4SJed Brown       pdd1 = dd;
26237045ce4SJed Brown       if (degrees[k] == j) {
26337045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
26437045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
26537045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
26637045ce4SJed Brown       }
26737045ce4SJed Brown     }
26837045ce4SJed Brown   }
26937045ce4SJed Brown   PetscFunctionReturn(0);
27037045ce4SJed Brown }
27137045ce4SJed Brown 
27237045ce4SJed Brown #undef __FUNCT__
27337045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature"
27437045ce4SJed Brown /*@
27537045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
27637045ce4SJed Brown 
27737045ce4SJed Brown    Not Collective
27837045ce4SJed Brown 
27937045ce4SJed Brown    Input Arguments:
28037045ce4SJed Brown +  npoints - number of points
28137045ce4SJed Brown .  a - left end of interval (often-1)
28237045ce4SJed Brown -  b - right end of interval (often +1)
28337045ce4SJed Brown 
28437045ce4SJed Brown    Output Arguments:
28537045ce4SJed Brown +  x - quadrature points
28637045ce4SJed Brown -  w - quadrature weights
28737045ce4SJed Brown 
28837045ce4SJed Brown    Level: intermediate
28937045ce4SJed Brown 
29037045ce4SJed Brown    References:
29137045ce4SJed Brown    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
29237045ce4SJed Brown 
29337045ce4SJed Brown .seealso: PetscDTLegendreEval()
29437045ce4SJed Brown @*/
29537045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
29637045ce4SJed Brown {
29737045ce4SJed Brown   PetscErrorCode ierr;
29837045ce4SJed Brown   PetscInt       i;
29937045ce4SJed Brown   PetscReal      *work;
30037045ce4SJed Brown   PetscScalar    *Z;
30137045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
30237045ce4SJed Brown 
30337045ce4SJed Brown   PetscFunctionBegin;
3040bfcf5a5SMatthew G. Knepley   ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr);
30537045ce4SJed Brown   /* Set up the Golub-Welsch system */
30637045ce4SJed Brown   for (i=0; i<npoints; i++) {
30737045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
30837045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
30937045ce4SJed Brown   }
310dcca6d9dSJed Brown   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
311c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
31237045ce4SJed Brown   LDZ  = N;
31337045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
3148b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
31537045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
3161c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
31737045ce4SJed Brown 
31837045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
31937045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
32037045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
32137045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
3220d644c17SKarl Rupp 
32388393a60SJed Brown     w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
32437045ce4SJed Brown   }
32537045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
32637045ce4SJed Brown   PetscFunctionReturn(0);
32737045ce4SJed Brown }
328194825f6SJed Brown 
329194825f6SJed Brown #undef __FUNCT__
330744bafbcSMatthew G. Knepley #define __FUNCT__ "PetscDTGaussTensorQuadrature"
331744bafbcSMatthew G. Knepley /*@
332744bafbcSMatthew G. Knepley   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
333744bafbcSMatthew G. Knepley 
334744bafbcSMatthew G. Knepley   Not Collective
335744bafbcSMatthew G. Knepley 
336744bafbcSMatthew G. Knepley   Input Arguments:
337744bafbcSMatthew G. Knepley + dim     - The spatial dimension
338744bafbcSMatthew G. Knepley . npoints - number of points in one dimension
339744bafbcSMatthew G. Knepley . a       - left end of interval (often-1)
340744bafbcSMatthew G. Knepley - b       - right end of interval (often +1)
341744bafbcSMatthew G. Knepley 
342744bafbcSMatthew G. Knepley   Output Argument:
343744bafbcSMatthew G. Knepley . q - A PetscQuadrature object
344744bafbcSMatthew G. Knepley 
345744bafbcSMatthew G. Knepley   Level: intermediate
346744bafbcSMatthew G. Knepley 
347744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
348744bafbcSMatthew G. Knepley @*/
349744bafbcSMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
350744bafbcSMatthew G. Knepley {
351744bafbcSMatthew G. Knepley   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k;
352744bafbcSMatthew G. Knepley   PetscReal     *x, *w, *xw, *ww;
353744bafbcSMatthew G. Knepley   PetscErrorCode ierr;
354744bafbcSMatthew G. Knepley 
355744bafbcSMatthew G. Knepley   PetscFunctionBegin;
356744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
357744bafbcSMatthew G. Knepley   ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr);
358744bafbcSMatthew G. Knepley   /* Set up the Golub-Welsch system */
359744bafbcSMatthew G. Knepley   switch (dim) {
360744bafbcSMatthew G. Knepley   case 0:
361744bafbcSMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
362744bafbcSMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
363744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
364744bafbcSMatthew G. Knepley     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
365744bafbcSMatthew G. Knepley     x[0] = 0.0;
366744bafbcSMatthew G. Knepley     w[0] = 1.0;
367744bafbcSMatthew G. Knepley     break;
368744bafbcSMatthew G. Knepley   case 1:
369744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr);
370744bafbcSMatthew G. Knepley     break;
371744bafbcSMatthew G. Knepley   case 2:
372744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
373744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
374744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
375744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
376744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+0] = xw[i];
377744bafbcSMatthew G. Knepley         x[(i*npoints+j)*dim+1] = xw[j];
378744bafbcSMatthew G. Knepley         w[i*npoints+j]         = ww[i] * ww[j];
379744bafbcSMatthew G. Knepley       }
380744bafbcSMatthew G. Knepley     }
381744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
382744bafbcSMatthew G. Knepley     break;
383744bafbcSMatthew G. Knepley   case 3:
384744bafbcSMatthew G. Knepley     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
385744bafbcSMatthew G. Knepley     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
386744bafbcSMatthew G. Knepley     for (i = 0; i < npoints; ++i) {
387744bafbcSMatthew G. Knepley       for (j = 0; j < npoints; ++j) {
388744bafbcSMatthew G. Knepley         for (k = 0; k < npoints; ++k) {
389744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
390744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
391744bafbcSMatthew G. Knepley           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
392744bafbcSMatthew G. Knepley           w[(i*npoints+j)*npoints+k]         = ww[i] * ww[j] * ww[k];
393744bafbcSMatthew G. Knepley         }
394744bafbcSMatthew G. Knepley       }
395744bafbcSMatthew G. Knepley     }
396744bafbcSMatthew G. Knepley     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
397744bafbcSMatthew G. Knepley     break;
398744bafbcSMatthew G. Knepley   default:
399744bafbcSMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
400744bafbcSMatthew G. Knepley   }
401744bafbcSMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
402744bafbcSMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr);
403744bafbcSMatthew G. Knepley   PetscFunctionReturn(0);
404744bafbcSMatthew G. Knepley }
405744bafbcSMatthew G. Knepley 
406744bafbcSMatthew G. Knepley #undef __FUNCT__
407494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal"
408494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
409494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
410494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
411494e7359SMatthew G. Knepley {
412494e7359SMatthew G. Knepley   PetscReal f = 1.0;
413494e7359SMatthew G. Knepley   PetscInt  i;
414494e7359SMatthew G. Knepley 
415494e7359SMatthew G. Knepley   PetscFunctionBegin;
416494e7359SMatthew G. Knepley   for (i = 1; i < n+1; ++i) f *= i;
417494e7359SMatthew G. Knepley   *factorial = f;
418494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
419494e7359SMatthew G. Knepley }
420494e7359SMatthew G. Knepley 
421494e7359SMatthew G. Knepley #undef __FUNCT__
422494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi"
423494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
424494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
425494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
426494e7359SMatthew G. Knepley {
427494e7359SMatthew G. Knepley   PetscReal apb, pn1, pn2;
428494e7359SMatthew G. Knepley   PetscInt  k;
429494e7359SMatthew G. Knepley 
430494e7359SMatthew G. Knepley   PetscFunctionBegin;
431494e7359SMatthew G. Knepley   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
432494e7359SMatthew G. Knepley   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
433494e7359SMatthew G. Knepley   apb = a + b;
434494e7359SMatthew G. Knepley   pn2 = 1.0;
435494e7359SMatthew G. Knepley   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
436494e7359SMatthew G. Knepley   *P  = 0.0;
437494e7359SMatthew G. Knepley   for (k = 2; k < n+1; ++k) {
438494e7359SMatthew G. Knepley     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
439494e7359SMatthew G. Knepley     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
440494e7359SMatthew G. Knepley     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
441494e7359SMatthew G. Knepley     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
442494e7359SMatthew G. Knepley 
443494e7359SMatthew G. Knepley     a2  = a2 / a1;
444494e7359SMatthew G. Knepley     a3  = a3 / a1;
445494e7359SMatthew G. Knepley     a4  = a4 / a1;
446494e7359SMatthew G. Knepley     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
447494e7359SMatthew G. Knepley     pn2 = pn1;
448494e7359SMatthew G. Knepley     pn1 = *P;
449494e7359SMatthew G. Knepley   }
450494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
451494e7359SMatthew G. Knepley }
452494e7359SMatthew G. Knepley 
453494e7359SMatthew G. Knepley #undef __FUNCT__
454494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative"
455494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
456494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
457494e7359SMatthew G. Knepley {
458494e7359SMatthew G. Knepley   PetscReal      nP;
459494e7359SMatthew G. Knepley   PetscErrorCode ierr;
460494e7359SMatthew G. Knepley 
461494e7359SMatthew G. Knepley   PetscFunctionBegin;
462494e7359SMatthew G. Knepley   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
463494e7359SMatthew G. Knepley   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
464494e7359SMatthew G. Knepley   *P   = 0.5 * (a + b + n + 1) * nP;
465494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
466494e7359SMatthew G. Knepley }
467494e7359SMatthew G. Knepley 
468494e7359SMatthew G. Knepley #undef __FUNCT__
469494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
470494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
471494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
472494e7359SMatthew G. Knepley {
473494e7359SMatthew G. Knepley   PetscFunctionBegin;
474494e7359SMatthew G. Knepley   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
475494e7359SMatthew G. Knepley   *eta = y;
476494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
477494e7359SMatthew G. Knepley }
478494e7359SMatthew G. Knepley 
479494e7359SMatthew G. Knepley #undef __FUNCT__
480494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
481494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
482494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
483494e7359SMatthew G. Knepley {
484494e7359SMatthew G. Knepley   PetscFunctionBegin;
485494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
486494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
487494e7359SMatthew G. Knepley   *zeta = z;
488494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
489494e7359SMatthew G. Knepley }
490494e7359SMatthew G. Knepley 
491494e7359SMatthew G. Knepley #undef __FUNCT__
492494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
493494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
494494e7359SMatthew G. Knepley {
495494e7359SMatthew G. Knepley   PetscInt       maxIter = 100;
496494e7359SMatthew G. Knepley   PetscReal      eps     = 1.0e-8;
497a8291ba1SSatish Balay   PetscReal      a1, a2, a3, a4, a5, a6;
498494e7359SMatthew G. Knepley   PetscInt       k;
499494e7359SMatthew G. Knepley   PetscErrorCode ierr;
500494e7359SMatthew G. Knepley 
501494e7359SMatthew G. Knepley   PetscFunctionBegin;
502a8291ba1SSatish Balay 
5038b49ba18SBarry Smith   a1      = PetscPowReal(2.0, a+b+1);
504a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA)
5050646a658SBarry Smith   a2      = PetscTGamma(a + npoints + 1);
5060646a658SBarry Smith   a3      = PetscTGamma(b + npoints + 1);
5070646a658SBarry Smith   a4      = PetscTGamma(a + b + npoints + 1);
508a8291ba1SSatish Balay #else
509a8291ba1SSatish Balay   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
510a8291ba1SSatish Balay #endif
511a8291ba1SSatish Balay 
512494e7359SMatthew G. Knepley   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
513494e7359SMatthew G. Knepley   a6   = a1 * a2 * a3 / a4 / a5;
514494e7359SMatthew G. Knepley   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
515494e7359SMatthew G. Knepley    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
516494e7359SMatthew G. Knepley   for (k = 0; k < npoints; ++k) {
5178b49ba18SBarry Smith     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
518494e7359SMatthew G. Knepley     PetscInt  j;
519494e7359SMatthew G. Knepley 
520494e7359SMatthew G. Knepley     if (k > 0) r = 0.5 * (r + x[k-1]);
521494e7359SMatthew G. Knepley     for (j = 0; j < maxIter; ++j) {
522494e7359SMatthew G. Knepley       PetscReal s = 0.0, delta, f, fp;
523494e7359SMatthew G. Knepley       PetscInt  i;
524494e7359SMatthew G. Knepley 
525494e7359SMatthew G. Knepley       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
526494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
527494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
528494e7359SMatthew G. Knepley       delta = f / (fp - f * s);
529494e7359SMatthew G. Knepley       r     = r - delta;
53077b4d14cSPeter Brune       if (PetscAbsReal(delta) < eps) break;
531494e7359SMatthew G. Knepley     }
532494e7359SMatthew G. Knepley     x[k] = r;
533494e7359SMatthew G. Knepley     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
534494e7359SMatthew G. Knepley     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
535494e7359SMatthew G. Knepley   }
536494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
537494e7359SMatthew G. Knepley }
538494e7359SMatthew G. Knepley 
539494e7359SMatthew G. Knepley #undef __FUNCT__
540494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
541fd9d31fbSMatthew G. Knepley /*@C
542494e7359SMatthew G. Knepley   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
543494e7359SMatthew G. Knepley 
544494e7359SMatthew G. Knepley   Not Collective
545494e7359SMatthew G. Knepley 
546494e7359SMatthew G. Knepley   Input Arguments:
547494e7359SMatthew G. Knepley + dim   - The simplex dimension
548744bafbcSMatthew G. Knepley . order - The number of points in one dimension
549494e7359SMatthew G. Knepley . a     - left end of interval (often-1)
550494e7359SMatthew G. Knepley - b     - right end of interval (often +1)
551494e7359SMatthew G. Knepley 
552744bafbcSMatthew G. Knepley   Output Argument:
553552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
554494e7359SMatthew G. Knepley 
555494e7359SMatthew G. Knepley   Level: intermediate
556494e7359SMatthew G. Knepley 
557494e7359SMatthew G. Knepley   References:
558494e7359SMatthew G. Knepley   Karniadakis and Sherwin.
559494e7359SMatthew G. Knepley   FIAT
560494e7359SMatthew G. Knepley 
561744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
562494e7359SMatthew G. Knepley @*/
563552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
564494e7359SMatthew G. Knepley {
565552aa4f7SMatthew G. Knepley   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
566494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
567494e7359SMatthew G. Knepley   PetscInt       i, j, k;
568494e7359SMatthew G. Knepley   PetscErrorCode ierr;
569494e7359SMatthew G. Knepley 
570494e7359SMatthew G. Knepley   PetscFunctionBegin;
571494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
572785e854fSJed Brown   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
573785e854fSJed Brown   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
574494e7359SMatthew G. Knepley   switch (dim) {
575707aa5c5SMatthew G. Knepley   case 0:
576707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
577707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
578785e854fSJed Brown     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
579785e854fSJed Brown     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
580707aa5c5SMatthew G. Knepley     x[0] = 0.0;
581707aa5c5SMatthew G. Knepley     w[0] = 1.0;
582707aa5c5SMatthew G. Knepley     break;
583494e7359SMatthew G. Knepley   case 1:
584552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
585494e7359SMatthew G. Knepley     break;
586494e7359SMatthew G. Knepley   case 2:
587dcca6d9dSJed Brown     ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr);
588552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
589552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
590552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
591552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
592552aa4f7SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
593552aa4f7SMatthew G. Knepley         w[i*order+j] = 0.5 * wx[i] * wy[j];
594494e7359SMatthew G. Knepley       }
595494e7359SMatthew G. Knepley     }
596494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
597494e7359SMatthew G. Knepley     break;
598494e7359SMatthew G. Knepley   case 3:
599dcca6d9dSJed Brown     ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr);
600552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
601552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
602552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
603552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
604552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
605552aa4f7SMatthew G. Knepley         for (k = 0; k < order; ++k) {
606552aa4f7SMatthew G. Knepley           ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*order+j)*order+k)*3+0], &x[((i*order+j)*order+k)*3+1], &x[((i*order+j)*order+k)*3+2]);CHKERRQ(ierr);
607552aa4f7SMatthew G. Knepley           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
608494e7359SMatthew G. Knepley         }
609494e7359SMatthew G. Knepley       }
610494e7359SMatthew G. Knepley     }
611494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
612494e7359SMatthew G. Knepley     break;
613494e7359SMatthew G. Knepley   default:
614494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
615494e7359SMatthew G. Knepley   }
61621454ff5SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
61721454ff5SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
618494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
619494e7359SMatthew G. Knepley }
620494e7359SMatthew G. Knepley 
621494e7359SMatthew G. Knepley #undef __FUNCT__
622194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR"
623194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
624194825f6SJed Brown  * A in column-major format
625194825f6SJed Brown  * Ainv in row-major format
626194825f6SJed Brown  * tau has length m
627194825f6SJed Brown  * worksize must be >= max(1,n)
628194825f6SJed Brown  */
629194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
630194825f6SJed Brown {
631194825f6SJed Brown   PetscErrorCode ierr;
632194825f6SJed Brown   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
633194825f6SJed Brown   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
634194825f6SJed Brown 
635194825f6SJed Brown   PetscFunctionBegin;
636194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
637194825f6SJed Brown   {
638194825f6SJed Brown     PetscInt i,j;
639dcca6d9dSJed Brown     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
640194825f6SJed Brown     for (j=0; j<n; j++) {
641194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
642194825f6SJed Brown     }
643194825f6SJed Brown     mstride = m;
644194825f6SJed Brown   }
645194825f6SJed Brown #else
646194825f6SJed Brown   A = A_in;
647194825f6SJed Brown   Ainv = Ainv_out;
648194825f6SJed Brown #endif
649194825f6SJed Brown 
650194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
651194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
652194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
653194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
654194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
655001a771dSBarry Smith   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
656194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
657194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
658194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
659194825f6SJed Brown 
660194825f6SJed Brown   /* Extract an explicit representation of Q */
661194825f6SJed Brown   Q = Ainv;
662194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
663194825f6SJed Brown   K = N;                        /* full rank */
664001a771dSBarry Smith   PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
665194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
666194825f6SJed Brown 
667194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
668194825f6SJed Brown   Alpha = 1.0;
669194825f6SJed Brown   ldb = lda;
670001a771dSBarry Smith   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
671194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
672194825f6SJed Brown 
673194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
674194825f6SJed Brown   {
675194825f6SJed Brown     PetscInt i;
676194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
677194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
678194825f6SJed Brown   }
679194825f6SJed Brown #endif
680194825f6SJed Brown   PetscFunctionReturn(0);
681194825f6SJed Brown }
682194825f6SJed Brown 
683194825f6SJed Brown #undef __FUNCT__
684194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate"
685194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
686194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
687194825f6SJed Brown {
688194825f6SJed Brown   PetscErrorCode ierr;
689194825f6SJed Brown   PetscReal      *Bv;
690194825f6SJed Brown   PetscInt       i,j;
691194825f6SJed Brown 
692194825f6SJed Brown   PetscFunctionBegin;
693785e854fSJed Brown   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
694194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
695194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
696194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
697194825f6SJed Brown   for (i=0; i<ninterval; i++) {
698194825f6SJed Brown     for (j=0; j<ndegree; j++) {
699194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
700194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
701194825f6SJed Brown     }
702194825f6SJed Brown   }
703194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
704194825f6SJed Brown   PetscFunctionReturn(0);
705194825f6SJed Brown }
706194825f6SJed Brown 
707194825f6SJed Brown #undef __FUNCT__
708194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly"
709194825f6SJed Brown /*@
710194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
711194825f6SJed Brown 
712194825f6SJed Brown    Not Collective
713194825f6SJed Brown 
714194825f6SJed Brown    Input Arguments:
715194825f6SJed Brown +  degree - degree of reconstruction polynomial
716194825f6SJed Brown .  nsource - number of source intervals
717194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
718194825f6SJed Brown .  ntarget - number of target intervals
719194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
720194825f6SJed Brown 
721194825f6SJed Brown    Output Arguments:
722194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
723194825f6SJed Brown 
724194825f6SJed Brown    Level: advanced
725194825f6SJed Brown 
726194825f6SJed Brown .seealso: PetscDTLegendreEval()
727194825f6SJed Brown @*/
728194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
729194825f6SJed Brown {
730194825f6SJed Brown   PetscErrorCode ierr;
731194825f6SJed Brown   PetscInt       i,j,k,*bdegrees,worksize;
732194825f6SJed Brown   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
733194825f6SJed Brown   PetscScalar    *tau,*work;
734194825f6SJed Brown 
735194825f6SJed Brown   PetscFunctionBegin;
736194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
737194825f6SJed Brown   PetscValidRealPointer(targetx,5);
738194825f6SJed Brown   PetscValidRealPointer(R,6);
739194825f6SJed 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);
740194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
741194825f6SJed Brown   for (i=0; i<nsource; i++) {
74257622a8eSBarry 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]);
743194825f6SJed Brown   }
744194825f6SJed Brown   for (i=0; i<ntarget; i++) {
74557622a8eSBarry 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]);
746194825f6SJed Brown   }
747194825f6SJed Brown #endif
748194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
749194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
750194825f6SJed Brown   center = (xmin + xmax)/2;
751194825f6SJed Brown   hscale = (xmax - xmin)/2;
752194825f6SJed Brown   worksize = nsource;
753dcca6d9dSJed Brown   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
754dcca6d9dSJed Brown   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
755194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
756194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
757194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
758194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
759194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
760194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
761194825f6SJed Brown   for (i=0; i<ntarget; i++) {
762194825f6SJed Brown     PetscReal rowsum = 0;
763194825f6SJed Brown     for (j=0; j<nsource; j++) {
764194825f6SJed Brown       PetscReal sum = 0;
765194825f6SJed Brown       for (k=0; k<degree+1; k++) {
766194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
767194825f6SJed Brown       }
768194825f6SJed Brown       R[i*nsource+j] = sum;
769194825f6SJed Brown       rowsum += sum;
770194825f6SJed Brown     }
771194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
772194825f6SJed Brown   }
773194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
774194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
775194825f6SJed Brown   PetscFunctionReturn(0);
776194825f6SJed Brown }
777