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> 11*21454ff5SMatthew 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__ 17*21454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureCreate" 18*21454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 19*21454ff5SMatthew G. Knepley { 20*21454ff5SMatthew G. Knepley PetscErrorCode ierr; 21*21454ff5SMatthew G. Knepley 22*21454ff5SMatthew G. Knepley PetscFunctionBegin; 23*21454ff5SMatthew G. Knepley PetscValidPointer(q, 2); 24*21454ff5SMatthew G. Knepley ierr = DMInitializePackage();CHKERRQ(ierr); 25*21454ff5SMatthew G. Knepley ierr = PetscHeaderCreate(*q,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 26*21454ff5SMatthew G. Knepley (*q)->dim = -1; 27*21454ff5SMatthew G. Knepley (*q)->numPoints = 0; 28*21454ff5SMatthew G. Knepley (*q)->points = NULL; 29*21454ff5SMatthew G. Knepley (*q)->weights = NULL; 30*21454ff5SMatthew G. Knepley PetscFunctionReturn(0); 31*21454ff5SMatthew G. Knepley } 32*21454ff5SMatthew G. Knepley 33*21454ff5SMatthew 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; 40*21454ff5SMatthew G. Knepley if (!*q) PetscFunctionReturn(0); 41*21454ff5SMatthew G. Knepley PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 42*21454ff5SMatthew G. Knepley if (--((PetscObject)(*q))->refct > 0) { 43*21454ff5SMatthew G. Knepley *q = NULL; 44*21454ff5SMatthew G. Knepley PetscFunctionReturn(0); 45*21454ff5SMatthew G. Knepley } 46*21454ff5SMatthew G. Knepley ierr = PetscFree((*q)->points);CHKERRQ(ierr); 47*21454ff5SMatthew G. Knepley ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 48*21454ff5SMatthew G. Knepley ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 49*21454ff5SMatthew G. Knepley PetscFunctionReturn(0); 50*21454ff5SMatthew G. Knepley } 51*21454ff5SMatthew G. Knepley 52*21454ff5SMatthew G. Knepley #undef __FUNCT__ 53*21454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetData" 54*21454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 55*21454ff5SMatthew G. Knepley { 56*21454ff5SMatthew G. Knepley PetscFunctionBegin; 57*21454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 58*21454ff5SMatthew G. Knepley if (dim) { 59*21454ff5SMatthew G. Knepley PetscValidPointer(dim, 2); 60*21454ff5SMatthew G. Knepley *dim = q->dim; 61*21454ff5SMatthew G. Knepley } 62*21454ff5SMatthew G. Knepley if (npoints) { 63*21454ff5SMatthew G. Knepley PetscValidPointer(npoints, 3); 64*21454ff5SMatthew G. Knepley *npoints = q->numPoints; 65*21454ff5SMatthew G. Knepley } 66*21454ff5SMatthew G. Knepley if (points) { 67*21454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 68*21454ff5SMatthew G. Knepley *points = q->points; 69*21454ff5SMatthew G. Knepley } 70*21454ff5SMatthew G. Knepley if (weights) { 71*21454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 72*21454ff5SMatthew G. Knepley *weights = q->weights; 73*21454ff5SMatthew G. Knepley } 74*21454ff5SMatthew G. Knepley PetscFunctionReturn(0); 75*21454ff5SMatthew G. Knepley } 76*21454ff5SMatthew G. Knepley 77*21454ff5SMatthew G. Knepley #undef __FUNCT__ 78*21454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetData" 79*21454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 80*21454ff5SMatthew G. Knepley { 81*21454ff5SMatthew G. Knepley PetscFunctionBegin; 82*21454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 83*21454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 84*21454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 85*21454ff5SMatthew G. Knepley if (points) { 86*21454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 87*21454ff5SMatthew G. Knepley q->points = points; 88*21454ff5SMatthew G. Knepley } 89*21454ff5SMatthew G. Knepley if (weights) { 90*21454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 91*21454ff5SMatthew G. Knepley q->weights = weights; 92*21454ff5SMatthew 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; 104*21454ff5SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n (", quad->numPoints);CHKERRQ(ierr); 105*21454ff5SMatthew G. Knepley for (q = 0; q < quad->numPoints; ++q) { 106*21454ff5SMatthew G. Knepley for (d = 0; d < quad->dim; ++d) { 107f9fd7fdbSMatthew G. Knepley if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 108*21454ff5SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%g\n", quad->points[q*quad->dim+d]);CHKERRQ(ierr); 109f9fd7fdbSMatthew G. Knepley } 110*21454ff5SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", 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 } 21937045ce4SJed Brown ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 220dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 221c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 22237045ce4SJed Brown LDZ = N; 22337045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 2248b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 22537045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 2261c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 22737045ce4SJed Brown 22837045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 22937045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 23037045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 23137045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 2320d644c17SKarl Rupp 23337045ce4SJed Brown w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints])); 23437045ce4SJed Brown } 23537045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 23637045ce4SJed Brown PetscFunctionReturn(0); 23737045ce4SJed Brown } 238194825f6SJed Brown 239194825f6SJed Brown #undef __FUNCT__ 240494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal" 241494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 242494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 243494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 244494e7359SMatthew G. Knepley { 245494e7359SMatthew G. Knepley PetscReal f = 1.0; 246494e7359SMatthew G. Knepley PetscInt i; 247494e7359SMatthew G. Knepley 248494e7359SMatthew G. Knepley PetscFunctionBegin; 249494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 250494e7359SMatthew G. Knepley *factorial = f; 251494e7359SMatthew G. Knepley PetscFunctionReturn(0); 252494e7359SMatthew G. Knepley } 253494e7359SMatthew G. Knepley 254494e7359SMatthew G. Knepley #undef __FUNCT__ 255494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi" 256494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 257494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 258494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 259494e7359SMatthew G. Knepley { 260494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 261494e7359SMatthew G. Knepley PetscInt k; 262494e7359SMatthew G. Knepley 263494e7359SMatthew G. Knepley PetscFunctionBegin; 264494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 265494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 266494e7359SMatthew G. Knepley apb = a + b; 267494e7359SMatthew G. Knepley pn2 = 1.0; 268494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 269494e7359SMatthew G. Knepley *P = 0.0; 270494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 271494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 272494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 273494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 274494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 275494e7359SMatthew G. Knepley 276494e7359SMatthew G. Knepley a2 = a2 / a1; 277494e7359SMatthew G. Knepley a3 = a3 / a1; 278494e7359SMatthew G. Knepley a4 = a4 / a1; 279494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 280494e7359SMatthew G. Knepley pn2 = pn1; 281494e7359SMatthew G. Knepley pn1 = *P; 282494e7359SMatthew G. Knepley } 283494e7359SMatthew G. Knepley PetscFunctionReturn(0); 284494e7359SMatthew G. Knepley } 285494e7359SMatthew G. Knepley 286494e7359SMatthew G. Knepley #undef __FUNCT__ 287494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative" 288494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 289494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 290494e7359SMatthew G. Knepley { 291494e7359SMatthew G. Knepley PetscReal nP; 292494e7359SMatthew G. Knepley PetscErrorCode ierr; 293494e7359SMatthew G. Knepley 294494e7359SMatthew G. Knepley PetscFunctionBegin; 295494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 296494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 297494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 298494e7359SMatthew G. Knepley PetscFunctionReturn(0); 299494e7359SMatthew G. Knepley } 300494e7359SMatthew G. Knepley 301494e7359SMatthew G. Knepley #undef __FUNCT__ 302494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 303494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 304494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 305494e7359SMatthew G. Knepley { 306494e7359SMatthew G. Knepley PetscFunctionBegin; 307494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 308494e7359SMatthew G. Knepley *eta = y; 309494e7359SMatthew G. Knepley PetscFunctionReturn(0); 310494e7359SMatthew G. Knepley } 311494e7359SMatthew G. Knepley 312494e7359SMatthew G. Knepley #undef __FUNCT__ 313494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 314494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 315494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 316494e7359SMatthew G. Knepley { 317494e7359SMatthew G. Knepley PetscFunctionBegin; 318494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 319494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 320494e7359SMatthew G. Knepley *zeta = z; 321494e7359SMatthew G. Knepley PetscFunctionReturn(0); 322494e7359SMatthew G. Knepley } 323494e7359SMatthew G. Knepley 324494e7359SMatthew G. Knepley #undef __FUNCT__ 325494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 326494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 327494e7359SMatthew G. Knepley { 328494e7359SMatthew G. Knepley PetscInt maxIter = 100; 329494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 330a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 331494e7359SMatthew G. Knepley PetscInt k; 332494e7359SMatthew G. Knepley PetscErrorCode ierr; 333494e7359SMatthew G. Knepley 334494e7359SMatthew G. Knepley PetscFunctionBegin; 335a8291ba1SSatish Balay 3368b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 337a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 338a8291ba1SSatish Balay a2 = tgamma(a + npoints + 1); 339a8291ba1SSatish Balay a3 = tgamma(b + npoints + 1); 340a8291ba1SSatish Balay a4 = tgamma(a + b + npoints + 1); 341a8291ba1SSatish Balay #else 342a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 343a8291ba1SSatish Balay #endif 344a8291ba1SSatish Balay 345494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 346494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 347494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 348494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 349494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 3508b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 351494e7359SMatthew G. Knepley PetscInt j; 352494e7359SMatthew G. Knepley 353494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 354494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 355494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 356494e7359SMatthew G. Knepley PetscInt i; 357494e7359SMatthew G. Knepley 358494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 359494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 360494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 361494e7359SMatthew G. Knepley delta = f / (fp - f * s); 362494e7359SMatthew G. Knepley r = r - delta; 363494e7359SMatthew G. Knepley if (fabs(delta) < eps) break; 364494e7359SMatthew G. Knepley } 365494e7359SMatthew G. Knepley x[k] = r; 366494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 367494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 368494e7359SMatthew G. Knepley } 369494e7359SMatthew G. Knepley PetscFunctionReturn(0); 370494e7359SMatthew G. Knepley } 371494e7359SMatthew G. Knepley 372494e7359SMatthew G. Knepley #undef __FUNCT__ 373494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 374fd9d31fbSMatthew G. Knepley /*@C 375494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 376494e7359SMatthew G. Knepley 377494e7359SMatthew G. Knepley Not Collective 378494e7359SMatthew G. Knepley 379494e7359SMatthew G. Knepley Input Arguments: 380494e7359SMatthew G. Knepley + dim - The simplex dimension 381552aa4f7SMatthew G. Knepley . order - The quadrature order 382494e7359SMatthew G. Knepley . a - left end of interval (often-1) 383494e7359SMatthew G. Knepley - b - right end of interval (often +1) 384494e7359SMatthew G. Knepley 385494e7359SMatthew G. Knepley Output Arguments: 386552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 387494e7359SMatthew G. Knepley 388494e7359SMatthew G. Knepley Level: intermediate 389494e7359SMatthew G. Knepley 390494e7359SMatthew G. Knepley References: 391494e7359SMatthew G. Knepley Karniadakis and Sherwin. 392494e7359SMatthew G. Knepley FIAT 393494e7359SMatthew G. Knepley 394494e7359SMatthew G. Knepley .seealso: PetscDTGaussQuadrature() 395494e7359SMatthew G. Knepley @*/ 396552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q) 397494e7359SMatthew G. Knepley { 398552aa4f7SMatthew G. Knepley PetscInt npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order; 399494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 400494e7359SMatthew G. Knepley PetscInt i, j, k; 401494e7359SMatthew G. Knepley PetscErrorCode ierr; 402494e7359SMatthew G. Knepley 403494e7359SMatthew G. Knepley PetscFunctionBegin; 404494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 405785e854fSJed Brown ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 406785e854fSJed Brown ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 407494e7359SMatthew G. Knepley switch (dim) { 408707aa5c5SMatthew G. Knepley case 0: 409707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 410707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 411785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 412785e854fSJed Brown ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 413707aa5c5SMatthew G. Knepley x[0] = 0.0; 414707aa5c5SMatthew G. Knepley w[0] = 1.0; 415707aa5c5SMatthew G. Knepley break; 416494e7359SMatthew G. Knepley case 1: 417552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr); 418494e7359SMatthew G. Knepley break; 419494e7359SMatthew G. Knepley case 2: 420dcca6d9dSJed Brown ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr); 421552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 422552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 423552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 424552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 425552aa4f7SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr); 426552aa4f7SMatthew G. Knepley w[i*order+j] = 0.5 * wx[i] * wy[j]; 427494e7359SMatthew G. Knepley } 428494e7359SMatthew G. Knepley } 429494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 430494e7359SMatthew G. Knepley break; 431494e7359SMatthew G. Knepley case 3: 432dcca6d9dSJed Brown ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr); 433552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 434552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 435552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 436552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 437552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 438552aa4f7SMatthew G. Knepley for (k = 0; k < order; ++k) { 439552aa4f7SMatthew 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); 440552aa4f7SMatthew G. Knepley w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k]; 441494e7359SMatthew G. Knepley } 442494e7359SMatthew G. Knepley } 443494e7359SMatthew G. Knepley } 444494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 445494e7359SMatthew G. Knepley break; 446494e7359SMatthew G. Knepley default: 447494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 448494e7359SMatthew G. Knepley } 449*21454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 450*21454ff5SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr); 451494e7359SMatthew G. Knepley PetscFunctionReturn(0); 452494e7359SMatthew G. Knepley } 453494e7359SMatthew G. Knepley 454494e7359SMatthew G. Knepley #undef __FUNCT__ 455194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR" 456194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 457194825f6SJed Brown * A in column-major format 458194825f6SJed Brown * Ainv in row-major format 459194825f6SJed Brown * tau has length m 460194825f6SJed Brown * worksize must be >= max(1,n) 461194825f6SJed Brown */ 462194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 463194825f6SJed Brown { 464194825f6SJed Brown PetscErrorCode ierr; 465194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 466194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 467194825f6SJed Brown 468194825f6SJed Brown PetscFunctionBegin; 469194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 470194825f6SJed Brown { 471194825f6SJed Brown PetscInt i,j; 472dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 473194825f6SJed Brown for (j=0; j<n; j++) { 474194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 475194825f6SJed Brown } 476194825f6SJed Brown mstride = m; 477194825f6SJed Brown } 478194825f6SJed Brown #else 479194825f6SJed Brown A = A_in; 480194825f6SJed Brown Ainv = Ainv_out; 481194825f6SJed Brown #endif 482194825f6SJed Brown 483194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 484194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 485194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 486194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 487194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 488194825f6SJed Brown LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info); 489194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 490194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 491194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 492194825f6SJed Brown 493194825f6SJed Brown /* Extract an explicit representation of Q */ 494194825f6SJed Brown Q = Ainv; 495194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 496194825f6SJed Brown K = N; /* full rank */ 497194825f6SJed Brown LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info); 498194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 499194825f6SJed Brown 500194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 501194825f6SJed Brown Alpha = 1.0; 502194825f6SJed Brown ldb = lda; 503194825f6SJed Brown BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb); 504194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 505194825f6SJed Brown 506194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 507194825f6SJed Brown { 508194825f6SJed Brown PetscInt i; 509194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 510194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 511194825f6SJed Brown } 512194825f6SJed Brown #endif 513194825f6SJed Brown PetscFunctionReturn(0); 514194825f6SJed Brown } 515194825f6SJed Brown 516194825f6SJed Brown #undef __FUNCT__ 517194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate" 518194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 519194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 520194825f6SJed Brown { 521194825f6SJed Brown PetscErrorCode ierr; 522194825f6SJed Brown PetscReal *Bv; 523194825f6SJed Brown PetscInt i,j; 524194825f6SJed Brown 525194825f6SJed Brown PetscFunctionBegin; 526785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 527194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 528194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 529194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 530194825f6SJed Brown for (i=0; i<ninterval; i++) { 531194825f6SJed Brown for (j=0; j<ndegree; j++) { 532194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 533194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 534194825f6SJed Brown } 535194825f6SJed Brown } 536194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 537194825f6SJed Brown PetscFunctionReturn(0); 538194825f6SJed Brown } 539194825f6SJed Brown 540194825f6SJed Brown #undef __FUNCT__ 541194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly" 542194825f6SJed Brown /*@ 543194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 544194825f6SJed Brown 545194825f6SJed Brown Not Collective 546194825f6SJed Brown 547194825f6SJed Brown Input Arguments: 548194825f6SJed Brown + degree - degree of reconstruction polynomial 549194825f6SJed Brown . nsource - number of source intervals 550194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 551194825f6SJed Brown . ntarget - number of target intervals 552194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 553194825f6SJed Brown 554194825f6SJed Brown Output Arguments: 555194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 556194825f6SJed Brown 557194825f6SJed Brown Level: advanced 558194825f6SJed Brown 559194825f6SJed Brown .seealso: PetscDTLegendreEval() 560194825f6SJed Brown @*/ 561194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 562194825f6SJed Brown { 563194825f6SJed Brown PetscErrorCode ierr; 564194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 565194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 566194825f6SJed Brown PetscScalar *tau,*work; 567194825f6SJed Brown 568194825f6SJed Brown PetscFunctionBegin; 569194825f6SJed Brown PetscValidRealPointer(sourcex,3); 570194825f6SJed Brown PetscValidRealPointer(targetx,5); 571194825f6SJed Brown PetscValidRealPointer(R,6); 572194825f6SJed 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); 573194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 574194825f6SJed Brown for (i=0; i<nsource; i++) { 57557622a8eSBarry 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]); 576194825f6SJed Brown } 577194825f6SJed Brown for (i=0; i<ntarget; i++) { 57857622a8eSBarry 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]); 579194825f6SJed Brown } 580194825f6SJed Brown #endif 581194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 582194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 583194825f6SJed Brown center = (xmin + xmax)/2; 584194825f6SJed Brown hscale = (xmax - xmin)/2; 585194825f6SJed Brown worksize = nsource; 586dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 587dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 588194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 589194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 590194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 591194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 592194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 593194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 594194825f6SJed Brown for (i=0; i<ntarget; i++) { 595194825f6SJed Brown PetscReal rowsum = 0; 596194825f6SJed Brown for (j=0; j<nsource; j++) { 597194825f6SJed Brown PetscReal sum = 0; 598194825f6SJed Brown for (k=0; k<degree+1; k++) { 599194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 600194825f6SJed Brown } 601194825f6SJed Brown R[i*nsource+j] = sum; 602194825f6SJed Brown rowsum += sum; 603194825f6SJed Brown } 604194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 605194825f6SJed Brown } 606194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 607194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 608194825f6SJed Brown PetscFunctionReturn(0); 609194825f6SJed Brown } 610