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" 2821454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 2921454ff5SMatthew G. Knepley { 3021454ff5SMatthew G. Knepley PetscErrorCode ierr; 3121454ff5SMatthew G. Knepley 3221454ff5SMatthew G. Knepley PetscFunctionBegin; 3321454ff5SMatthew G. Knepley PetscValidPointer(q, 2); 3421454ff5SMatthew G. Knepley ierr = DMInitializePackage();CHKERRQ(ierr); 3521454ff5SMatthew G. Knepley ierr = PetscHeaderCreate(*q,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 3621454ff5SMatthew G. Knepley (*q)->dim = -1; 3721454ff5SMatthew G. Knepley (*q)->numPoints = 0; 3821454ff5SMatthew G. Knepley (*q)->points = NULL; 3921454ff5SMatthew G. Knepley (*q)->weights = NULL; 4021454ff5SMatthew G. Knepley PetscFunctionReturn(0); 4121454ff5SMatthew G. Knepley } 4221454ff5SMatthew G. Knepley 4321454ff5SMatthew G. Knepley #undef __FUNCT__ 44bfa639d9SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDestroy" 45bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 46bfa639d9SMatthew G. Knepley { 47bfa639d9SMatthew G. Knepley PetscErrorCode ierr; 48bfa639d9SMatthew G. Knepley 49bfa639d9SMatthew G. Knepley PetscFunctionBegin; 5021454ff5SMatthew G. Knepley if (!*q) PetscFunctionReturn(0); 5121454ff5SMatthew G. Knepley PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 5221454ff5SMatthew G. Knepley if (--((PetscObject)(*q))->refct > 0) { 5321454ff5SMatthew G. Knepley *q = NULL; 5421454ff5SMatthew G. Knepley PetscFunctionReturn(0); 5521454ff5SMatthew G. Knepley } 5621454ff5SMatthew G. Knepley ierr = PetscFree((*q)->points);CHKERRQ(ierr); 5721454ff5SMatthew G. Knepley ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 5821454ff5SMatthew G. Knepley ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 5921454ff5SMatthew G. Knepley PetscFunctionReturn(0); 6021454ff5SMatthew G. Knepley } 6121454ff5SMatthew G. Knepley 6221454ff5SMatthew G. Knepley #undef __FUNCT__ 6321454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetData" 6421454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 6521454ff5SMatthew G. Knepley { 6621454ff5SMatthew G. Knepley PetscFunctionBegin; 6721454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 6821454ff5SMatthew G. Knepley if (dim) { 6921454ff5SMatthew G. Knepley PetscValidPointer(dim, 2); 7021454ff5SMatthew G. Knepley *dim = q->dim; 7121454ff5SMatthew G. Knepley } 7221454ff5SMatthew G. Knepley if (npoints) { 7321454ff5SMatthew G. Knepley PetscValidPointer(npoints, 3); 7421454ff5SMatthew G. Knepley *npoints = q->numPoints; 7521454ff5SMatthew G. Knepley } 7621454ff5SMatthew G. Knepley if (points) { 7721454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 7821454ff5SMatthew G. Knepley *points = q->points; 7921454ff5SMatthew G. Knepley } 8021454ff5SMatthew G. Knepley if (weights) { 8121454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 8221454ff5SMatthew G. Knepley *weights = q->weights; 8321454ff5SMatthew G. Knepley } 8421454ff5SMatthew G. Knepley PetscFunctionReturn(0); 8521454ff5SMatthew G. Knepley } 8621454ff5SMatthew G. Knepley 8721454ff5SMatthew G. Knepley #undef __FUNCT__ 8821454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetData" 8921454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 9021454ff5SMatthew G. Knepley { 9121454ff5SMatthew G. Knepley PetscFunctionBegin; 9221454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 9321454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 9421454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 9521454ff5SMatthew G. Knepley if (points) { 9621454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 9721454ff5SMatthew G. Knepley q->points = points; 9821454ff5SMatthew G. Knepley } 9921454ff5SMatthew G. Knepley if (weights) { 10021454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 10121454ff5SMatthew G. Knepley q->weights = weights; 10221454ff5SMatthew G. Knepley } 103f9fd7fdbSMatthew G. Knepley PetscFunctionReturn(0); 104f9fd7fdbSMatthew G. Knepley } 105f9fd7fdbSMatthew G. Knepley 106f9fd7fdbSMatthew G. Knepley #undef __FUNCT__ 107f9fd7fdbSMatthew G. Knepley #define __FUNCT__ "PetscQuadratureView" 108f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 109f9fd7fdbSMatthew G. Knepley { 110f9fd7fdbSMatthew G. Knepley PetscInt q, d; 111f9fd7fdbSMatthew G. Knepley PetscErrorCode ierr; 112f9fd7fdbSMatthew G. Knepley 113f9fd7fdbSMatthew G. Knepley PetscFunctionBegin; 114*98c3331eSBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr); 11521454ff5SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n (", quad->numPoints);CHKERRQ(ierr); 11621454ff5SMatthew G. Knepley for (q = 0; q < quad->numPoints; ++q) { 11721454ff5SMatthew G. Knepley for (d = 0; d < quad->dim; ++d) { 118f9fd7fdbSMatthew G. Knepley if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 119ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 120f9fd7fdbSMatthew G. Knepley } 121ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr); 122f9fd7fdbSMatthew G. Knepley } 123bfa639d9SMatthew G. Knepley PetscFunctionReturn(0); 124bfa639d9SMatthew G. Knepley } 125bfa639d9SMatthew G. Knepley 126bfa639d9SMatthew G. Knepley #undef __FUNCT__ 12737045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval" 12837045ce4SJed Brown /*@ 12937045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 13037045ce4SJed Brown 13137045ce4SJed Brown Not Collective 13237045ce4SJed Brown 13337045ce4SJed Brown Input Arguments: 13437045ce4SJed Brown + npoints - number of spatial points to evaluate at 13537045ce4SJed Brown . points - array of locations to evaluate at 13637045ce4SJed Brown . ndegree - number of basis degrees to evaluate 13737045ce4SJed Brown - degrees - sorted array of degrees to evaluate 13837045ce4SJed Brown 13937045ce4SJed Brown Output Arguments: 1400298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 1410298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 1420298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 14337045ce4SJed Brown 14437045ce4SJed Brown Level: intermediate 14537045ce4SJed Brown 14637045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 14737045ce4SJed Brown @*/ 14837045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 14937045ce4SJed Brown { 15037045ce4SJed Brown PetscInt i,maxdegree; 15137045ce4SJed Brown 15237045ce4SJed Brown PetscFunctionBegin; 15337045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 15437045ce4SJed Brown maxdegree = degrees[ndegree-1]; 15537045ce4SJed Brown for (i=0; i<npoints; i++) { 15637045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 15737045ce4SJed Brown PetscInt j,k; 15837045ce4SJed Brown x = points[i]; 15937045ce4SJed Brown pm2 = 0; 16037045ce4SJed Brown pm1 = 1; 16137045ce4SJed Brown pd2 = 0; 16237045ce4SJed Brown pd1 = 0; 16337045ce4SJed Brown pdd2 = 0; 16437045ce4SJed Brown pdd1 = 0; 16537045ce4SJed Brown k = 0; 16637045ce4SJed Brown if (degrees[k] == 0) { 16737045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 16837045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 16937045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 17037045ce4SJed Brown k++; 17137045ce4SJed Brown } 17237045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 17337045ce4SJed Brown PetscReal p,d,dd; 17437045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 17537045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 17637045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 17737045ce4SJed Brown pm2 = pm1; 17837045ce4SJed Brown pm1 = p; 17937045ce4SJed Brown pd2 = pd1; 18037045ce4SJed Brown pd1 = d; 18137045ce4SJed Brown pdd2 = pdd1; 18237045ce4SJed Brown pdd1 = dd; 18337045ce4SJed Brown if (degrees[k] == j) { 18437045ce4SJed Brown if (B) B[i*ndegree+k] = p; 18537045ce4SJed Brown if (D) D[i*ndegree+k] = d; 18637045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 18737045ce4SJed Brown } 18837045ce4SJed Brown } 18937045ce4SJed Brown } 19037045ce4SJed Brown PetscFunctionReturn(0); 19137045ce4SJed Brown } 19237045ce4SJed Brown 19337045ce4SJed Brown #undef __FUNCT__ 19437045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature" 19537045ce4SJed Brown /*@ 19637045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 19737045ce4SJed Brown 19837045ce4SJed Brown Not Collective 19937045ce4SJed Brown 20037045ce4SJed Brown Input Arguments: 20137045ce4SJed Brown + npoints - number of points 20237045ce4SJed Brown . a - left end of interval (often-1) 20337045ce4SJed Brown - b - right end of interval (often +1) 20437045ce4SJed Brown 20537045ce4SJed Brown Output Arguments: 20637045ce4SJed Brown + x - quadrature points 20737045ce4SJed Brown - w - quadrature weights 20837045ce4SJed Brown 20937045ce4SJed Brown Level: intermediate 21037045ce4SJed Brown 21137045ce4SJed Brown References: 21237045ce4SJed Brown Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 21337045ce4SJed Brown 21437045ce4SJed Brown .seealso: PetscDTLegendreEval() 21537045ce4SJed Brown @*/ 21637045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 21737045ce4SJed Brown { 21837045ce4SJed Brown PetscErrorCode ierr; 21937045ce4SJed Brown PetscInt i; 22037045ce4SJed Brown PetscReal *work; 22137045ce4SJed Brown PetscScalar *Z; 22237045ce4SJed Brown PetscBLASInt N,LDZ,info; 22337045ce4SJed Brown 22437045ce4SJed Brown PetscFunctionBegin; 2250bfcf5a5SMatthew G. Knepley ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 22637045ce4SJed Brown /* Set up the Golub-Welsch system */ 22737045ce4SJed Brown for (i=0; i<npoints; i++) { 22837045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 22937045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 23037045ce4SJed Brown } 231dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 232c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 23337045ce4SJed Brown LDZ = N; 23437045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 2358b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 23637045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 2371c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 23837045ce4SJed Brown 23937045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 24037045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 24137045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 24237045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 2430d644c17SKarl Rupp 24488393a60SJed Brown w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 24537045ce4SJed Brown } 24637045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 24737045ce4SJed Brown PetscFunctionReturn(0); 24837045ce4SJed Brown } 249194825f6SJed Brown 250194825f6SJed Brown #undef __FUNCT__ 251744bafbcSMatthew G. Knepley #define __FUNCT__ "PetscDTGaussTensorQuadrature" 252744bafbcSMatthew G. Knepley /*@ 253744bafbcSMatthew G. Knepley PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 254744bafbcSMatthew G. Knepley 255744bafbcSMatthew G. Knepley Not Collective 256744bafbcSMatthew G. Knepley 257744bafbcSMatthew G. Knepley Input Arguments: 258744bafbcSMatthew G. Knepley + dim - The spatial dimension 259744bafbcSMatthew G. Knepley . npoints - number of points in one dimension 260744bafbcSMatthew G. Knepley . a - left end of interval (often-1) 261744bafbcSMatthew G. Knepley - b - right end of interval (often +1) 262744bafbcSMatthew G. Knepley 263744bafbcSMatthew G. Knepley Output Argument: 264744bafbcSMatthew G. Knepley . q - A PetscQuadrature object 265744bafbcSMatthew G. Knepley 266744bafbcSMatthew G. Knepley Level: intermediate 267744bafbcSMatthew G. Knepley 268744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 269744bafbcSMatthew G. Knepley @*/ 270744bafbcSMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 271744bafbcSMatthew G. Knepley { 272744bafbcSMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k; 273744bafbcSMatthew G. Knepley PetscReal *x, *w, *xw, *ww; 274744bafbcSMatthew G. Knepley PetscErrorCode ierr; 275744bafbcSMatthew G. Knepley 276744bafbcSMatthew G. Knepley PetscFunctionBegin; 277744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 278744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr); 279744bafbcSMatthew G. Knepley /* Set up the Golub-Welsch system */ 280744bafbcSMatthew G. Knepley switch (dim) { 281744bafbcSMatthew G. Knepley case 0: 282744bafbcSMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 283744bafbcSMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 284744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 285744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 286744bafbcSMatthew G. Knepley x[0] = 0.0; 287744bafbcSMatthew G. Knepley w[0] = 1.0; 288744bafbcSMatthew G. Knepley break; 289744bafbcSMatthew G. Knepley case 1: 290744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr); 291744bafbcSMatthew G. Knepley break; 292744bafbcSMatthew G. Knepley case 2: 293744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 294744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 295744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 296744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 297744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+0] = xw[i]; 298744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+1] = xw[j]; 299744bafbcSMatthew G. Knepley w[i*npoints+j] = ww[i] * ww[j]; 300744bafbcSMatthew G. Knepley } 301744bafbcSMatthew G. Knepley } 302744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 303744bafbcSMatthew G. Knepley break; 304744bafbcSMatthew G. Knepley case 3: 305744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 306744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 307744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 308744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 309744bafbcSMatthew G. Knepley for (k = 0; k < npoints; ++k) { 310744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 311744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 312744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 313744bafbcSMatthew G. Knepley w[(i*npoints+j)*npoints+k] = ww[i] * ww[j] * ww[k]; 314744bafbcSMatthew G. Knepley } 315744bafbcSMatthew G. Knepley } 316744bafbcSMatthew G. Knepley } 317744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 318744bafbcSMatthew G. Knepley break; 319744bafbcSMatthew G. Knepley default: 320744bafbcSMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 321744bafbcSMatthew G. Knepley } 322744bafbcSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 323744bafbcSMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr); 324744bafbcSMatthew G. Knepley PetscFunctionReturn(0); 325744bafbcSMatthew G. Knepley } 326744bafbcSMatthew G. Knepley 327744bafbcSMatthew G. Knepley #undef __FUNCT__ 328494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal" 329494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 330494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 331494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 332494e7359SMatthew G. Knepley { 333494e7359SMatthew G. Knepley PetscReal f = 1.0; 334494e7359SMatthew G. Knepley PetscInt i; 335494e7359SMatthew G. Knepley 336494e7359SMatthew G. Knepley PetscFunctionBegin; 337494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 338494e7359SMatthew G. Knepley *factorial = f; 339494e7359SMatthew G. Knepley PetscFunctionReturn(0); 340494e7359SMatthew G. Knepley } 341494e7359SMatthew G. Knepley 342494e7359SMatthew G. Knepley #undef __FUNCT__ 343494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi" 344494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 345494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 346494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 347494e7359SMatthew G. Knepley { 348494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 349494e7359SMatthew G. Knepley PetscInt k; 350494e7359SMatthew G. Knepley 351494e7359SMatthew G. Knepley PetscFunctionBegin; 352494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 353494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 354494e7359SMatthew G. Knepley apb = a + b; 355494e7359SMatthew G. Knepley pn2 = 1.0; 356494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 357494e7359SMatthew G. Knepley *P = 0.0; 358494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 359494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 360494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 361494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 362494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 363494e7359SMatthew G. Knepley 364494e7359SMatthew G. Knepley a2 = a2 / a1; 365494e7359SMatthew G. Knepley a3 = a3 / a1; 366494e7359SMatthew G. Knepley a4 = a4 / a1; 367494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 368494e7359SMatthew G. Knepley pn2 = pn1; 369494e7359SMatthew G. Knepley pn1 = *P; 370494e7359SMatthew G. Knepley } 371494e7359SMatthew G. Knepley PetscFunctionReturn(0); 372494e7359SMatthew G. Knepley } 373494e7359SMatthew G. Knepley 374494e7359SMatthew G. Knepley #undef __FUNCT__ 375494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative" 376494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 377494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 378494e7359SMatthew G. Knepley { 379494e7359SMatthew G. Knepley PetscReal nP; 380494e7359SMatthew G. Knepley PetscErrorCode ierr; 381494e7359SMatthew G. Knepley 382494e7359SMatthew G. Knepley PetscFunctionBegin; 383494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 384494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 385494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 386494e7359SMatthew G. Knepley PetscFunctionReturn(0); 387494e7359SMatthew G. Knepley } 388494e7359SMatthew G. Knepley 389494e7359SMatthew G. Knepley #undef __FUNCT__ 390494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 391494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 392494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 393494e7359SMatthew G. Knepley { 394494e7359SMatthew G. Knepley PetscFunctionBegin; 395494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 396494e7359SMatthew G. Knepley *eta = y; 397494e7359SMatthew G. Knepley PetscFunctionReturn(0); 398494e7359SMatthew G. Knepley } 399494e7359SMatthew G. Knepley 400494e7359SMatthew G. Knepley #undef __FUNCT__ 401494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 402494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 403494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 404494e7359SMatthew G. Knepley { 405494e7359SMatthew G. Knepley PetscFunctionBegin; 406494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 407494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 408494e7359SMatthew G. Knepley *zeta = z; 409494e7359SMatthew G. Knepley PetscFunctionReturn(0); 410494e7359SMatthew G. Knepley } 411494e7359SMatthew G. Knepley 412494e7359SMatthew G. Knepley #undef __FUNCT__ 413494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 414494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 415494e7359SMatthew G. Knepley { 416494e7359SMatthew G. Knepley PetscInt maxIter = 100; 417494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 418a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 419494e7359SMatthew G. Knepley PetscInt k; 420494e7359SMatthew G. Knepley PetscErrorCode ierr; 421494e7359SMatthew G. Knepley 422494e7359SMatthew G. Knepley PetscFunctionBegin; 423a8291ba1SSatish Balay 4248b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 425a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 4260646a658SBarry Smith a2 = PetscTGamma(a + npoints + 1); 4270646a658SBarry Smith a3 = PetscTGamma(b + npoints + 1); 4280646a658SBarry Smith a4 = PetscTGamma(a + b + npoints + 1); 429a8291ba1SSatish Balay #else 430a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 431a8291ba1SSatish Balay #endif 432a8291ba1SSatish Balay 433494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 434494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 435494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 436494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 437494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 4388b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 439494e7359SMatthew G. Knepley PetscInt j; 440494e7359SMatthew G. Knepley 441494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 442494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 443494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 444494e7359SMatthew G. Knepley PetscInt i; 445494e7359SMatthew G. Knepley 446494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 447494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 448494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 449494e7359SMatthew G. Knepley delta = f / (fp - f * s); 450494e7359SMatthew G. Knepley r = r - delta; 451001a771dSBarry Smith if (PetscAbs(delta) < eps) break; 452494e7359SMatthew G. Knepley } 453494e7359SMatthew G. Knepley x[k] = r; 454494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 455494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 456494e7359SMatthew G. Knepley } 457494e7359SMatthew G. Knepley PetscFunctionReturn(0); 458494e7359SMatthew G. Knepley } 459494e7359SMatthew G. Knepley 460494e7359SMatthew G. Knepley #undef __FUNCT__ 461494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 462fd9d31fbSMatthew G. Knepley /*@C 463494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 464494e7359SMatthew G. Knepley 465494e7359SMatthew G. Knepley Not Collective 466494e7359SMatthew G. Knepley 467494e7359SMatthew G. Knepley Input Arguments: 468494e7359SMatthew G. Knepley + dim - The simplex dimension 469744bafbcSMatthew G. Knepley . order - The number of points in one dimension 470494e7359SMatthew G. Knepley . a - left end of interval (often-1) 471494e7359SMatthew G. Knepley - b - right end of interval (often +1) 472494e7359SMatthew G. Knepley 473744bafbcSMatthew G. Knepley Output Argument: 474552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 475494e7359SMatthew G. Knepley 476494e7359SMatthew G. Knepley Level: intermediate 477494e7359SMatthew G. Knepley 478494e7359SMatthew G. Knepley References: 479494e7359SMatthew G. Knepley Karniadakis and Sherwin. 480494e7359SMatthew G. Knepley FIAT 481494e7359SMatthew G. Knepley 482744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 483494e7359SMatthew G. Knepley @*/ 484552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q) 485494e7359SMatthew G. Knepley { 486552aa4f7SMatthew G. Knepley PetscInt npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order; 487494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 488494e7359SMatthew G. Knepley PetscInt i, j, k; 489494e7359SMatthew G. Knepley PetscErrorCode ierr; 490494e7359SMatthew G. Knepley 491494e7359SMatthew G. Knepley PetscFunctionBegin; 492494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 493785e854fSJed Brown ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 494785e854fSJed Brown ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 495494e7359SMatthew G. Knepley switch (dim) { 496707aa5c5SMatthew G. Knepley case 0: 497707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 498707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 499785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 500785e854fSJed Brown ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 501707aa5c5SMatthew G. Knepley x[0] = 0.0; 502707aa5c5SMatthew G. Knepley w[0] = 1.0; 503707aa5c5SMatthew G. Knepley break; 504494e7359SMatthew G. Knepley case 1: 505552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr); 506494e7359SMatthew G. Knepley break; 507494e7359SMatthew G. Knepley case 2: 508dcca6d9dSJed Brown ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr); 509552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 510552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 511552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 512552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 513552aa4f7SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr); 514552aa4f7SMatthew G. Knepley w[i*order+j] = 0.5 * wx[i] * wy[j]; 515494e7359SMatthew G. Knepley } 516494e7359SMatthew G. Knepley } 517494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 518494e7359SMatthew G. Knepley break; 519494e7359SMatthew G. Knepley case 3: 520dcca6d9dSJed Brown ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr); 521552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 522552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 523552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 524552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 525552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 526552aa4f7SMatthew G. Knepley for (k = 0; k < order; ++k) { 527552aa4f7SMatthew 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); 528552aa4f7SMatthew G. Knepley w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k]; 529494e7359SMatthew G. Knepley } 530494e7359SMatthew G. Knepley } 531494e7359SMatthew G. Knepley } 532494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 533494e7359SMatthew G. Knepley break; 534494e7359SMatthew G. Knepley default: 535494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 536494e7359SMatthew G. Knepley } 53721454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 53821454ff5SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr); 539494e7359SMatthew G. Knepley PetscFunctionReturn(0); 540494e7359SMatthew G. Knepley } 541494e7359SMatthew G. Knepley 542494e7359SMatthew G. Knepley #undef __FUNCT__ 543194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR" 544194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 545194825f6SJed Brown * A in column-major format 546194825f6SJed Brown * Ainv in row-major format 547194825f6SJed Brown * tau has length m 548194825f6SJed Brown * worksize must be >= max(1,n) 549194825f6SJed Brown */ 550194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 551194825f6SJed Brown { 552194825f6SJed Brown PetscErrorCode ierr; 553194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 554194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 555194825f6SJed Brown 556194825f6SJed Brown PetscFunctionBegin; 557194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 558194825f6SJed Brown { 559194825f6SJed Brown PetscInt i,j; 560dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 561194825f6SJed Brown for (j=0; j<n; j++) { 562194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 563194825f6SJed Brown } 564194825f6SJed Brown mstride = m; 565194825f6SJed Brown } 566194825f6SJed Brown #else 567194825f6SJed Brown A = A_in; 568194825f6SJed Brown Ainv = Ainv_out; 569194825f6SJed Brown #endif 570194825f6SJed Brown 571194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 572194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 573194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 574194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 575194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 576001a771dSBarry Smith PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 577194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 578194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 579194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 580194825f6SJed Brown 581194825f6SJed Brown /* Extract an explicit representation of Q */ 582194825f6SJed Brown Q = Ainv; 583194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 584194825f6SJed Brown K = N; /* full rank */ 585001a771dSBarry Smith PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 586194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 587194825f6SJed Brown 588194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 589194825f6SJed Brown Alpha = 1.0; 590194825f6SJed Brown ldb = lda; 591001a771dSBarry Smith PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 592194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 593194825f6SJed Brown 594194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 595194825f6SJed Brown { 596194825f6SJed Brown PetscInt i; 597194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 598194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 599194825f6SJed Brown } 600194825f6SJed Brown #endif 601194825f6SJed Brown PetscFunctionReturn(0); 602194825f6SJed Brown } 603194825f6SJed Brown 604194825f6SJed Brown #undef __FUNCT__ 605194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate" 606194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 607194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 608194825f6SJed Brown { 609194825f6SJed Brown PetscErrorCode ierr; 610194825f6SJed Brown PetscReal *Bv; 611194825f6SJed Brown PetscInt i,j; 612194825f6SJed Brown 613194825f6SJed Brown PetscFunctionBegin; 614785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 615194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 616194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 617194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 618194825f6SJed Brown for (i=0; i<ninterval; i++) { 619194825f6SJed Brown for (j=0; j<ndegree; j++) { 620194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 621194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 622194825f6SJed Brown } 623194825f6SJed Brown } 624194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 625194825f6SJed Brown PetscFunctionReturn(0); 626194825f6SJed Brown } 627194825f6SJed Brown 628194825f6SJed Brown #undef __FUNCT__ 629194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly" 630194825f6SJed Brown /*@ 631194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 632194825f6SJed Brown 633194825f6SJed Brown Not Collective 634194825f6SJed Brown 635194825f6SJed Brown Input Arguments: 636194825f6SJed Brown + degree - degree of reconstruction polynomial 637194825f6SJed Brown . nsource - number of source intervals 638194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 639194825f6SJed Brown . ntarget - number of target intervals 640194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 641194825f6SJed Brown 642194825f6SJed Brown Output Arguments: 643194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 644194825f6SJed Brown 645194825f6SJed Brown Level: advanced 646194825f6SJed Brown 647194825f6SJed Brown .seealso: PetscDTLegendreEval() 648194825f6SJed Brown @*/ 649194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 650194825f6SJed Brown { 651194825f6SJed Brown PetscErrorCode ierr; 652194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 653194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 654194825f6SJed Brown PetscScalar *tau,*work; 655194825f6SJed Brown 656194825f6SJed Brown PetscFunctionBegin; 657194825f6SJed Brown PetscValidRealPointer(sourcex,3); 658194825f6SJed Brown PetscValidRealPointer(targetx,5); 659194825f6SJed Brown PetscValidRealPointer(R,6); 660194825f6SJed 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); 661194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 662194825f6SJed Brown for (i=0; i<nsource; i++) { 66357622a8eSBarry 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]); 664194825f6SJed Brown } 665194825f6SJed Brown for (i=0; i<ntarget; i++) { 66657622a8eSBarry 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]); 667194825f6SJed Brown } 668194825f6SJed Brown #endif 669194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 670194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 671194825f6SJed Brown center = (xmin + xmax)/2; 672194825f6SJed Brown hscale = (xmax - xmin)/2; 673194825f6SJed Brown worksize = nsource; 674dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 675dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 676194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 677194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 678194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 679194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 680194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 681194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 682194825f6SJed Brown for (i=0; i<ntarget; i++) { 683194825f6SJed Brown PetscReal rowsum = 0; 684194825f6SJed Brown for (j=0; j<nsource; j++) { 685194825f6SJed Brown PetscReal sum = 0; 686194825f6SJed Brown for (k=0; k<degree+1; k++) { 687194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 688194825f6SJed Brown } 689194825f6SJed Brown R[i*nsource+j] = sum; 690194825f6SJed Brown rowsum += sum; 691194825f6SJed Brown } 692194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 693194825f6SJed Brown } 694194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 695194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 696194825f6SJed Brown PetscFunctionReturn(0); 697194825f6SJed Brown } 698