xref: /petsc/src/dm/dt/interface/dt.c (revision ab15ae4336b482d5290728009e85a2f7ef5673ff)
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 
1637045ce4SJed Brown #undef __FUNCT__
1721454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureCreate"
1821454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
1921454ff5SMatthew G. Knepley {
2021454ff5SMatthew G. Knepley   PetscErrorCode ierr;
2121454ff5SMatthew G. Knepley 
2221454ff5SMatthew G. Knepley   PetscFunctionBegin;
2321454ff5SMatthew G. Knepley   PetscValidPointer(q, 2);
2421454ff5SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
2521454ff5SMatthew G. Knepley   ierr = PetscHeaderCreate(*q,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
2621454ff5SMatthew G. Knepley   (*q)->dim       = -1;
2721454ff5SMatthew G. Knepley   (*q)->numPoints = 0;
2821454ff5SMatthew G. Knepley   (*q)->points    = NULL;
2921454ff5SMatthew G. Knepley   (*q)->weights   = NULL;
3021454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
3121454ff5SMatthew G. Knepley }
3221454ff5SMatthew G. Knepley 
3321454ff5SMatthew G. Knepley #undef __FUNCT__
34bfa639d9SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDestroy"
35bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
36bfa639d9SMatthew G. Knepley {
37bfa639d9SMatthew G. Knepley   PetscErrorCode ierr;
38bfa639d9SMatthew G. Knepley 
39bfa639d9SMatthew G. Knepley   PetscFunctionBegin;
4021454ff5SMatthew G. Knepley   if (!*q) PetscFunctionReturn(0);
4121454ff5SMatthew G. Knepley   PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1);
4221454ff5SMatthew G. Knepley   if (--((PetscObject)(*q))->refct > 0) {
4321454ff5SMatthew G. Knepley     *q = NULL;
4421454ff5SMatthew G. Knepley     PetscFunctionReturn(0);
4521454ff5SMatthew G. Knepley   }
4621454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
4721454ff5SMatthew G. Knepley   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
4821454ff5SMatthew G. Knepley   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
4921454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
5021454ff5SMatthew G. Knepley }
5121454ff5SMatthew G. Knepley 
5221454ff5SMatthew G. Knepley #undef __FUNCT__
5321454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetData"
5421454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
5521454ff5SMatthew G. Knepley {
5621454ff5SMatthew G. Knepley   PetscFunctionBegin;
5721454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
5821454ff5SMatthew G. Knepley   if (dim) {
5921454ff5SMatthew G. Knepley     PetscValidPointer(dim, 2);
6021454ff5SMatthew G. Knepley     *dim = q->dim;
6121454ff5SMatthew G. Knepley   }
6221454ff5SMatthew G. Knepley   if (npoints) {
6321454ff5SMatthew G. Knepley     PetscValidPointer(npoints, 3);
6421454ff5SMatthew G. Knepley     *npoints = q->numPoints;
6521454ff5SMatthew G. Knepley   }
6621454ff5SMatthew G. Knepley   if (points) {
6721454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
6821454ff5SMatthew G. Knepley     *points = q->points;
6921454ff5SMatthew G. Knepley   }
7021454ff5SMatthew G. Knepley   if (weights) {
7121454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
7221454ff5SMatthew G. Knepley     *weights = q->weights;
7321454ff5SMatthew G. Knepley   }
7421454ff5SMatthew G. Knepley   PetscFunctionReturn(0);
7521454ff5SMatthew G. Knepley }
7621454ff5SMatthew G. Knepley 
7721454ff5SMatthew G. Knepley #undef __FUNCT__
7821454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetData"
7921454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
8021454ff5SMatthew G. Knepley {
8121454ff5SMatthew G. Knepley   PetscFunctionBegin;
8221454ff5SMatthew G. Knepley   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
8321454ff5SMatthew G. Knepley   if (dim >= 0)     q->dim       = dim;
8421454ff5SMatthew G. Knepley   if (npoints >= 0) q->numPoints = npoints;
8521454ff5SMatthew G. Knepley   if (points) {
8621454ff5SMatthew G. Knepley     PetscValidPointer(points, 4);
8721454ff5SMatthew G. Knepley     q->points = points;
8821454ff5SMatthew G. Knepley   }
8921454ff5SMatthew G. Knepley   if (weights) {
9021454ff5SMatthew G. Knepley     PetscValidPointer(weights, 5);
9121454ff5SMatthew G. Knepley     q->weights = weights;
9221454ff5SMatthew G. Knepley   }
93f9fd7fdbSMatthew G. Knepley   PetscFunctionReturn(0);
94f9fd7fdbSMatthew G. Knepley }
95f9fd7fdbSMatthew G. Knepley 
96f9fd7fdbSMatthew G. Knepley #undef __FUNCT__
97f9fd7fdbSMatthew G. Knepley #define __FUNCT__ "PetscQuadratureView"
98f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
99f9fd7fdbSMatthew G. Knepley {
100f9fd7fdbSMatthew G. Knepley   PetscInt       q, d;
101f9fd7fdbSMatthew G. Knepley   PetscErrorCode ierr;
102f9fd7fdbSMatthew G. Knepley 
103f9fd7fdbSMatthew G. Knepley   PetscFunctionBegin;
10421454ff5SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n  (", quad->numPoints);CHKERRQ(ierr);
10521454ff5SMatthew G. Knepley   for (q = 0; q < quad->numPoints; ++q) {
10621454ff5SMatthew G. Knepley     for (d = 0; d < quad->dim; ++d) {
107f9fd7fdbSMatthew G. Knepley       if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);
108*ab15ae43SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
109f9fd7fdbSMatthew G. Knepley     }
110*ab15ae43SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr);
111f9fd7fdbSMatthew G. Knepley   }
112bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
113bfa639d9SMatthew G. Knepley }
114bfa639d9SMatthew G. Knepley 
115bfa639d9SMatthew G. Knepley #undef __FUNCT__
11637045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval"
11737045ce4SJed Brown /*@
11837045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
11937045ce4SJed Brown 
12037045ce4SJed Brown    Not Collective
12137045ce4SJed Brown 
12237045ce4SJed Brown    Input Arguments:
12337045ce4SJed Brown +  npoints - number of spatial points to evaluate at
12437045ce4SJed Brown .  points - array of locations to evaluate at
12537045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
12637045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
12737045ce4SJed Brown 
12837045ce4SJed Brown    Output Arguments:
1290298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
1300298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
1310298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
13237045ce4SJed Brown 
13337045ce4SJed Brown    Level: intermediate
13437045ce4SJed Brown 
13537045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
13637045ce4SJed Brown @*/
13737045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
13837045ce4SJed Brown {
13937045ce4SJed Brown   PetscInt i,maxdegree;
14037045ce4SJed Brown 
14137045ce4SJed Brown   PetscFunctionBegin;
14237045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
14337045ce4SJed Brown   maxdegree = degrees[ndegree-1];
14437045ce4SJed Brown   for (i=0; i<npoints; i++) {
14537045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
14637045ce4SJed Brown     PetscInt  j,k;
14737045ce4SJed Brown     x    = points[i];
14837045ce4SJed Brown     pm2  = 0;
14937045ce4SJed Brown     pm1  = 1;
15037045ce4SJed Brown     pd2  = 0;
15137045ce4SJed Brown     pd1  = 0;
15237045ce4SJed Brown     pdd2 = 0;
15337045ce4SJed Brown     pdd1 = 0;
15437045ce4SJed Brown     k    = 0;
15537045ce4SJed Brown     if (degrees[k] == 0) {
15637045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
15737045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
15837045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
15937045ce4SJed Brown       k++;
16037045ce4SJed Brown     }
16137045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
16237045ce4SJed Brown       PetscReal p,d,dd;
16337045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
16437045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
16537045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
16637045ce4SJed Brown       pm2  = pm1;
16737045ce4SJed Brown       pm1  = p;
16837045ce4SJed Brown       pd2  = pd1;
16937045ce4SJed Brown       pd1  = d;
17037045ce4SJed Brown       pdd2 = pdd1;
17137045ce4SJed Brown       pdd1 = dd;
17237045ce4SJed Brown       if (degrees[k] == j) {
17337045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
17437045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
17537045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
17637045ce4SJed Brown       }
17737045ce4SJed Brown     }
17837045ce4SJed Brown   }
17937045ce4SJed Brown   PetscFunctionReturn(0);
18037045ce4SJed Brown }
18137045ce4SJed Brown 
18237045ce4SJed Brown #undef __FUNCT__
18337045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature"
18437045ce4SJed Brown /*@
18537045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
18637045ce4SJed Brown 
18737045ce4SJed Brown    Not Collective
18837045ce4SJed Brown 
18937045ce4SJed Brown    Input Arguments:
19037045ce4SJed Brown +  npoints - number of points
19137045ce4SJed Brown .  a - left end of interval (often-1)
19237045ce4SJed Brown -  b - right end of interval (often +1)
19337045ce4SJed Brown 
19437045ce4SJed Brown    Output Arguments:
19537045ce4SJed Brown +  x - quadrature points
19637045ce4SJed Brown -  w - quadrature weights
19737045ce4SJed Brown 
19837045ce4SJed Brown    Level: intermediate
19937045ce4SJed Brown 
20037045ce4SJed Brown    References:
20137045ce4SJed Brown    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
20237045ce4SJed Brown 
20337045ce4SJed Brown .seealso: PetscDTLegendreEval()
20437045ce4SJed Brown @*/
20537045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
20637045ce4SJed Brown {
20737045ce4SJed Brown   PetscErrorCode ierr;
20837045ce4SJed Brown   PetscInt       i;
20937045ce4SJed Brown   PetscReal      *work;
21037045ce4SJed Brown   PetscScalar    *Z;
21137045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
21237045ce4SJed Brown 
21337045ce4SJed Brown   PetscFunctionBegin;
21437045ce4SJed Brown   /* Set up the Golub-Welsch system */
21537045ce4SJed Brown   for (i=0; i<npoints; i++) {
21637045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
21737045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
21837045ce4SJed Brown   }
219dcca6d9dSJed Brown   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
220c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
22137045ce4SJed Brown   LDZ  = N;
22237045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
2238b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
22437045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
2251c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
22637045ce4SJed Brown 
22737045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
22837045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
22937045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
23037045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
2310d644c17SKarl Rupp 
23237045ce4SJed Brown     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
23337045ce4SJed Brown   }
23437045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
23537045ce4SJed Brown   PetscFunctionReturn(0);
23637045ce4SJed Brown }
237194825f6SJed Brown 
238194825f6SJed Brown #undef __FUNCT__
239494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal"
240494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
241494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
242494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
243494e7359SMatthew G. Knepley {
244494e7359SMatthew G. Knepley   PetscReal f = 1.0;
245494e7359SMatthew G. Knepley   PetscInt  i;
246494e7359SMatthew G. Knepley 
247494e7359SMatthew G. Knepley   PetscFunctionBegin;
248494e7359SMatthew G. Knepley   for (i = 1; i < n+1; ++i) f *= i;
249494e7359SMatthew G. Knepley   *factorial = f;
250494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
251494e7359SMatthew G. Knepley }
252494e7359SMatthew G. Knepley 
253494e7359SMatthew G. Knepley #undef __FUNCT__
254494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi"
255494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
256494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
257494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
258494e7359SMatthew G. Knepley {
259494e7359SMatthew G. Knepley   PetscReal apb, pn1, pn2;
260494e7359SMatthew G. Knepley   PetscInt  k;
261494e7359SMatthew G. Knepley 
262494e7359SMatthew G. Knepley   PetscFunctionBegin;
263494e7359SMatthew G. Knepley   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
264494e7359SMatthew G. Knepley   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
265494e7359SMatthew G. Knepley   apb = a + b;
266494e7359SMatthew G. Knepley   pn2 = 1.0;
267494e7359SMatthew G. Knepley   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
268494e7359SMatthew G. Knepley   *P  = 0.0;
269494e7359SMatthew G. Knepley   for (k = 2; k < n+1; ++k) {
270494e7359SMatthew G. Knepley     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
271494e7359SMatthew G. Knepley     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
272494e7359SMatthew G. Knepley     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
273494e7359SMatthew G. Knepley     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
274494e7359SMatthew G. Knepley 
275494e7359SMatthew G. Knepley     a2  = a2 / a1;
276494e7359SMatthew G. Knepley     a3  = a3 / a1;
277494e7359SMatthew G. Knepley     a4  = a4 / a1;
278494e7359SMatthew G. Knepley     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
279494e7359SMatthew G. Knepley     pn2 = pn1;
280494e7359SMatthew G. Knepley     pn1 = *P;
281494e7359SMatthew G. Knepley   }
282494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
283494e7359SMatthew G. Knepley }
284494e7359SMatthew G. Knepley 
285494e7359SMatthew G. Knepley #undef __FUNCT__
286494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative"
287494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
288494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
289494e7359SMatthew G. Knepley {
290494e7359SMatthew G. Knepley   PetscReal      nP;
291494e7359SMatthew G. Knepley   PetscErrorCode ierr;
292494e7359SMatthew G. Knepley 
293494e7359SMatthew G. Knepley   PetscFunctionBegin;
294494e7359SMatthew G. Knepley   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
295494e7359SMatthew G. Knepley   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
296494e7359SMatthew G. Knepley   *P   = 0.5 * (a + b + n + 1) * nP;
297494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
298494e7359SMatthew G. Knepley }
299494e7359SMatthew G. Knepley 
300494e7359SMatthew G. Knepley #undef __FUNCT__
301494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
302494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
303494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
304494e7359SMatthew G. Knepley {
305494e7359SMatthew G. Knepley   PetscFunctionBegin;
306494e7359SMatthew G. Knepley   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
307494e7359SMatthew G. Knepley   *eta = y;
308494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
309494e7359SMatthew G. Knepley }
310494e7359SMatthew G. Knepley 
311494e7359SMatthew G. Knepley #undef __FUNCT__
312494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
313494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
314494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
315494e7359SMatthew G. Knepley {
316494e7359SMatthew G. Knepley   PetscFunctionBegin;
317494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
318494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
319494e7359SMatthew G. Knepley   *zeta = z;
320494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
321494e7359SMatthew G. Knepley }
322494e7359SMatthew G. Knepley 
323494e7359SMatthew G. Knepley #undef __FUNCT__
324494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
325494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
326494e7359SMatthew G. Knepley {
327494e7359SMatthew G. Knepley   PetscInt       maxIter = 100;
328494e7359SMatthew G. Knepley   PetscReal      eps     = 1.0e-8;
329a8291ba1SSatish Balay   PetscReal      a1, a2, a3, a4, a5, a6;
330494e7359SMatthew G. Knepley   PetscInt       k;
331494e7359SMatthew G. Knepley   PetscErrorCode ierr;
332494e7359SMatthew G. Knepley 
333494e7359SMatthew G. Knepley   PetscFunctionBegin;
334a8291ba1SSatish Balay 
3358b49ba18SBarry Smith   a1      = PetscPowReal(2.0, a+b+1);
336a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA)
3370646a658SBarry Smith   a2      = PetscTGamma(a + npoints + 1);
3380646a658SBarry Smith   a3      = PetscTGamma(b + npoints + 1);
3390646a658SBarry Smith   a4      = PetscTGamma(a + b + npoints + 1);
340a8291ba1SSatish Balay #else
341a8291ba1SSatish Balay   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
342a8291ba1SSatish Balay #endif
343a8291ba1SSatish Balay 
344494e7359SMatthew G. Knepley   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
345494e7359SMatthew G. Knepley   a6   = a1 * a2 * a3 / a4 / a5;
346494e7359SMatthew G. Knepley   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
347494e7359SMatthew G. Knepley    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
348494e7359SMatthew G. Knepley   for (k = 0; k < npoints; ++k) {
3498b49ba18SBarry Smith     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
350494e7359SMatthew G. Knepley     PetscInt  j;
351494e7359SMatthew G. Knepley 
352494e7359SMatthew G. Knepley     if (k > 0) r = 0.5 * (r + x[k-1]);
353494e7359SMatthew G. Knepley     for (j = 0; j < maxIter; ++j) {
354494e7359SMatthew G. Knepley       PetscReal s = 0.0, delta, f, fp;
355494e7359SMatthew G. Knepley       PetscInt  i;
356494e7359SMatthew G. Knepley 
357494e7359SMatthew G. Knepley       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
358494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
359494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
360494e7359SMatthew G. Knepley       delta = f / (fp - f * s);
361494e7359SMatthew G. Knepley       r     = r - delta;
362001a771dSBarry Smith       if (PetscAbs(delta) < eps) break;
363494e7359SMatthew G. Knepley     }
364494e7359SMatthew G. Knepley     x[k] = r;
365494e7359SMatthew G. Knepley     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
366494e7359SMatthew G. Knepley     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
367494e7359SMatthew G. Knepley   }
368494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
369494e7359SMatthew G. Knepley }
370494e7359SMatthew G. Knepley 
371494e7359SMatthew G. Knepley #undef __FUNCT__
372494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
373fd9d31fbSMatthew G. Knepley /*@C
374494e7359SMatthew G. Knepley   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
375494e7359SMatthew G. Knepley 
376494e7359SMatthew G. Knepley   Not Collective
377494e7359SMatthew G. Knepley 
378494e7359SMatthew G. Knepley   Input Arguments:
379494e7359SMatthew G. Knepley + dim - The simplex dimension
380552aa4f7SMatthew G. Knepley . order - The quadrature order
381494e7359SMatthew G. Knepley . a - left end of interval (often-1)
382494e7359SMatthew G. Knepley - b - right end of interval (often +1)
383494e7359SMatthew G. Knepley 
384494e7359SMatthew G. Knepley   Output Arguments:
385552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
386494e7359SMatthew G. Knepley 
387494e7359SMatthew G. Knepley   Level: intermediate
388494e7359SMatthew G. Knepley 
389494e7359SMatthew G. Knepley   References:
390494e7359SMatthew G. Knepley   Karniadakis and Sherwin.
391494e7359SMatthew G. Knepley   FIAT
392494e7359SMatthew G. Knepley 
393494e7359SMatthew G. Knepley .seealso: PetscDTGaussQuadrature()
394494e7359SMatthew G. Knepley @*/
395552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
396494e7359SMatthew G. Knepley {
397552aa4f7SMatthew G. Knepley   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
398494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
399494e7359SMatthew G. Knepley   PetscInt       i, j, k;
400494e7359SMatthew G. Knepley   PetscErrorCode ierr;
401494e7359SMatthew G. Knepley 
402494e7359SMatthew G. Knepley   PetscFunctionBegin;
403494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
404785e854fSJed Brown   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
405785e854fSJed Brown   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
406494e7359SMatthew G. Knepley   switch (dim) {
407707aa5c5SMatthew G. Knepley   case 0:
408707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
409707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
410785e854fSJed Brown     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
411785e854fSJed Brown     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
412707aa5c5SMatthew G. Knepley     x[0] = 0.0;
413707aa5c5SMatthew G. Knepley     w[0] = 1.0;
414707aa5c5SMatthew G. Knepley     break;
415494e7359SMatthew G. Knepley   case 1:
416552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
417494e7359SMatthew G. Knepley     break;
418494e7359SMatthew G. Knepley   case 2:
419dcca6d9dSJed Brown     ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr);
420552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
421552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
422552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
423552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
424552aa4f7SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
425552aa4f7SMatthew G. Knepley         w[i*order+j] = 0.5 * wx[i] * wy[j];
426494e7359SMatthew G. Knepley       }
427494e7359SMatthew G. Knepley     }
428494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
429494e7359SMatthew G. Knepley     break;
430494e7359SMatthew G. Knepley   case 3:
431dcca6d9dSJed Brown     ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr);
432552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
433552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
434552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
435552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
436552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
437552aa4f7SMatthew G. Knepley         for (k = 0; k < order; ++k) {
438552aa4f7SMatthew 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);
439552aa4f7SMatthew G. Knepley           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
440494e7359SMatthew G. Knepley         }
441494e7359SMatthew G. Knepley       }
442494e7359SMatthew G. Knepley     }
443494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
444494e7359SMatthew G. Knepley     break;
445494e7359SMatthew G. Knepley   default:
446494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
447494e7359SMatthew G. Knepley   }
44821454ff5SMatthew G. Knepley   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
44921454ff5SMatthew G. Knepley   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
450494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
451494e7359SMatthew G. Knepley }
452494e7359SMatthew G. Knepley 
453494e7359SMatthew G. Knepley #undef __FUNCT__
454194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR"
455194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
456194825f6SJed Brown  * A in column-major format
457194825f6SJed Brown  * Ainv in row-major format
458194825f6SJed Brown  * tau has length m
459194825f6SJed Brown  * worksize must be >= max(1,n)
460194825f6SJed Brown  */
461194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
462194825f6SJed Brown {
463194825f6SJed Brown   PetscErrorCode ierr;
464194825f6SJed Brown   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
465194825f6SJed Brown   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
466194825f6SJed Brown 
467194825f6SJed Brown   PetscFunctionBegin;
468194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
469194825f6SJed Brown   {
470194825f6SJed Brown     PetscInt i,j;
471dcca6d9dSJed Brown     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
472194825f6SJed Brown     for (j=0; j<n; j++) {
473194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
474194825f6SJed Brown     }
475194825f6SJed Brown     mstride = m;
476194825f6SJed Brown   }
477194825f6SJed Brown #else
478194825f6SJed Brown   A = A_in;
479194825f6SJed Brown   Ainv = Ainv_out;
480194825f6SJed Brown #endif
481194825f6SJed Brown 
482194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
483194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
484194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
485194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
486194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
487001a771dSBarry Smith   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
488194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
489194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
490194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
491194825f6SJed Brown 
492194825f6SJed Brown   /* Extract an explicit representation of Q */
493194825f6SJed Brown   Q = Ainv;
494194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
495194825f6SJed Brown   K = N;                        /* full rank */
496001a771dSBarry Smith   PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
497194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
498194825f6SJed Brown 
499194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
500194825f6SJed Brown   Alpha = 1.0;
501194825f6SJed Brown   ldb = lda;
502001a771dSBarry Smith   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
503194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
504194825f6SJed Brown 
505194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
506194825f6SJed Brown   {
507194825f6SJed Brown     PetscInt i;
508194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
509194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
510194825f6SJed Brown   }
511194825f6SJed Brown #endif
512194825f6SJed Brown   PetscFunctionReturn(0);
513194825f6SJed Brown }
514194825f6SJed Brown 
515194825f6SJed Brown #undef __FUNCT__
516194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate"
517194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
518194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
519194825f6SJed Brown {
520194825f6SJed Brown   PetscErrorCode ierr;
521194825f6SJed Brown   PetscReal      *Bv;
522194825f6SJed Brown   PetscInt       i,j;
523194825f6SJed Brown 
524194825f6SJed Brown   PetscFunctionBegin;
525785e854fSJed Brown   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
526194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
527194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
528194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
529194825f6SJed Brown   for (i=0; i<ninterval; i++) {
530194825f6SJed Brown     for (j=0; j<ndegree; j++) {
531194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
532194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
533194825f6SJed Brown     }
534194825f6SJed Brown   }
535194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
536194825f6SJed Brown   PetscFunctionReturn(0);
537194825f6SJed Brown }
538194825f6SJed Brown 
539194825f6SJed Brown #undef __FUNCT__
540194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly"
541194825f6SJed Brown /*@
542194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
543194825f6SJed Brown 
544194825f6SJed Brown    Not Collective
545194825f6SJed Brown 
546194825f6SJed Brown    Input Arguments:
547194825f6SJed Brown +  degree - degree of reconstruction polynomial
548194825f6SJed Brown .  nsource - number of source intervals
549194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
550194825f6SJed Brown .  ntarget - number of target intervals
551194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
552194825f6SJed Brown 
553194825f6SJed Brown    Output Arguments:
554194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
555194825f6SJed Brown 
556194825f6SJed Brown    Level: advanced
557194825f6SJed Brown 
558194825f6SJed Brown .seealso: PetscDTLegendreEval()
559194825f6SJed Brown @*/
560194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
561194825f6SJed Brown {
562194825f6SJed Brown   PetscErrorCode ierr;
563194825f6SJed Brown   PetscInt       i,j,k,*bdegrees,worksize;
564194825f6SJed Brown   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
565194825f6SJed Brown   PetscScalar    *tau,*work;
566194825f6SJed Brown 
567194825f6SJed Brown   PetscFunctionBegin;
568194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
569194825f6SJed Brown   PetscValidRealPointer(targetx,5);
570194825f6SJed Brown   PetscValidRealPointer(R,6);
571194825f6SJed 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);
572194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
573194825f6SJed Brown   for (i=0; i<nsource; i++) {
57457622a8eSBarry 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]);
575194825f6SJed Brown   }
576194825f6SJed Brown   for (i=0; i<ntarget; i++) {
57757622a8eSBarry 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]);
578194825f6SJed Brown   }
579194825f6SJed Brown #endif
580194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
581194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
582194825f6SJed Brown   center = (xmin + xmax)/2;
583194825f6SJed Brown   hscale = (xmax - xmin)/2;
584194825f6SJed Brown   worksize = nsource;
585dcca6d9dSJed Brown   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
586dcca6d9dSJed Brown   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
587194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
588194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
589194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
590194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
591194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
592194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
593194825f6SJed Brown   for (i=0; i<ntarget; i++) {
594194825f6SJed Brown     PetscReal rowsum = 0;
595194825f6SJed Brown     for (j=0; j<nsource; j++) {
596194825f6SJed Brown       PetscReal sum = 0;
597194825f6SJed Brown       for (k=0; k<degree+1; k++) {
598194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
599194825f6SJed Brown       }
600194825f6SJed Brown       R[i*nsource+j] = sum;
601194825f6SJed Brown       rowsum += sum;
602194825f6SJed Brown     }
603194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
604194825f6SJed Brown   }
605194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
606194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
607194825f6SJed Brown   PetscFunctionReturn(0);
608194825f6SJed Brown }
609