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; 37*bcede257SMatthew G. Knepley (*q)->order = -1; 3821454ff5SMatthew G. Knepley (*q)->numPoints = 0; 3921454ff5SMatthew G. Knepley (*q)->points = NULL; 4021454ff5SMatthew G. Knepley (*q)->weights = NULL; 4121454ff5SMatthew G. Knepley PetscFunctionReturn(0); 4221454ff5SMatthew G. Knepley } 4321454ff5SMatthew G. Knepley 4421454ff5SMatthew G. Knepley #undef __FUNCT__ 45bfa639d9SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDestroy" 46bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 47bfa639d9SMatthew G. Knepley { 48bfa639d9SMatthew G. Knepley PetscErrorCode ierr; 49bfa639d9SMatthew G. Knepley 50bfa639d9SMatthew G. Knepley PetscFunctionBegin; 5121454ff5SMatthew G. Knepley if (!*q) PetscFunctionReturn(0); 5221454ff5SMatthew G. Knepley PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 5321454ff5SMatthew G. Knepley if (--((PetscObject)(*q))->refct > 0) { 5421454ff5SMatthew G. Knepley *q = NULL; 5521454ff5SMatthew G. Knepley PetscFunctionReturn(0); 5621454ff5SMatthew G. Knepley } 5721454ff5SMatthew G. Knepley ierr = PetscFree((*q)->points);CHKERRQ(ierr); 5821454ff5SMatthew G. Knepley ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 5921454ff5SMatthew G. Knepley ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 6021454ff5SMatthew G. Knepley PetscFunctionReturn(0); 6121454ff5SMatthew G. Knepley } 6221454ff5SMatthew G. Knepley 6321454ff5SMatthew G. Knepley #undef __FUNCT__ 64*bcede257SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetOrder" 65*bcede257SMatthew G. Knepley /*@ 66*bcede257SMatthew G. Knepley PetscQuadratureGetOrder - Return the quadrature information 67*bcede257SMatthew G. Knepley 68*bcede257SMatthew G. Knepley Not collective 69*bcede257SMatthew G. Knepley 70*bcede257SMatthew G. Knepley Input Parameter: 71*bcede257SMatthew G. Knepley . q - The PetscQuadrature object 72*bcede257SMatthew G. Knepley 73*bcede257SMatthew G. Knepley Output Parameter: 74*bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 75*bcede257SMatthew G. Knepley 76*bcede257SMatthew G. Knepley Output Parameter: 77*bcede257SMatthew G. Knepley 78*bcede257SMatthew G. Knepley Level: intermediate 79*bcede257SMatthew G. Knepley 80*bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 81*bcede257SMatthew G. Knepley @*/ 82*bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 83*bcede257SMatthew G. Knepley { 84*bcede257SMatthew G. Knepley PetscFunctionBegin; 85*bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 86*bcede257SMatthew G. Knepley PetscValidPointer(order, 2); 87*bcede257SMatthew G. Knepley *order = q->order; 88*bcede257SMatthew G. Knepley PetscFunctionReturn(0); 89*bcede257SMatthew G. Knepley } 90*bcede257SMatthew G. Knepley 91*bcede257SMatthew G. Knepley #undef __FUNCT__ 92*bcede257SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetOrder" 93*bcede257SMatthew G. Knepley /*@ 94*bcede257SMatthew G. Knepley PetscQuadratureSetOrder - Return the quadrature information 95*bcede257SMatthew G. Knepley 96*bcede257SMatthew G. Knepley Not collective 97*bcede257SMatthew G. Knepley 98*bcede257SMatthew G. Knepley Input Parameters: 99*bcede257SMatthew G. Knepley + q - The PetscQuadrature object 100*bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 101*bcede257SMatthew G. Knepley 102*bcede257SMatthew G. Knepley Level: intermediate 103*bcede257SMatthew G. Knepley 104*bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 105*bcede257SMatthew G. Knepley @*/ 106*bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 107*bcede257SMatthew G. Knepley { 108*bcede257SMatthew G. Knepley PetscFunctionBegin; 109*bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 110*bcede257SMatthew G. Knepley q->order = order; 111*bcede257SMatthew G. Knepley PetscFunctionReturn(0); 112*bcede257SMatthew G. Knepley } 113*bcede257SMatthew G. Knepley 114*bcede257SMatthew G. Knepley #undef __FUNCT__ 11521454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetData" 11621454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 11721454ff5SMatthew G. Knepley { 11821454ff5SMatthew G. Knepley PetscFunctionBegin; 11921454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 12021454ff5SMatthew G. Knepley if (dim) { 12121454ff5SMatthew G. Knepley PetscValidPointer(dim, 2); 12221454ff5SMatthew G. Knepley *dim = q->dim; 12321454ff5SMatthew G. Knepley } 12421454ff5SMatthew G. Knepley if (npoints) { 12521454ff5SMatthew G. Knepley PetscValidPointer(npoints, 3); 12621454ff5SMatthew G. Knepley *npoints = q->numPoints; 12721454ff5SMatthew G. Knepley } 12821454ff5SMatthew G. Knepley if (points) { 12921454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 13021454ff5SMatthew G. Knepley *points = q->points; 13121454ff5SMatthew G. Knepley } 13221454ff5SMatthew G. Knepley if (weights) { 13321454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 13421454ff5SMatthew G. Knepley *weights = q->weights; 13521454ff5SMatthew G. Knepley } 13621454ff5SMatthew G. Knepley PetscFunctionReturn(0); 13721454ff5SMatthew G. Knepley } 13821454ff5SMatthew G. Knepley 13921454ff5SMatthew G. Knepley #undef __FUNCT__ 14021454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetData" 14121454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 14221454ff5SMatthew G. Knepley { 14321454ff5SMatthew G. Knepley PetscFunctionBegin; 14421454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 14521454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 14621454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 14721454ff5SMatthew G. Knepley if (points) { 14821454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 14921454ff5SMatthew G. Knepley q->points = points; 15021454ff5SMatthew G. Knepley } 15121454ff5SMatthew G. Knepley if (weights) { 15221454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 15321454ff5SMatthew G. Knepley q->weights = weights; 15421454ff5SMatthew G. Knepley } 155f9fd7fdbSMatthew G. Knepley PetscFunctionReturn(0); 156f9fd7fdbSMatthew G. Knepley } 157f9fd7fdbSMatthew G. Knepley 158f9fd7fdbSMatthew G. Knepley #undef __FUNCT__ 159f9fd7fdbSMatthew G. Knepley #define __FUNCT__ "PetscQuadratureView" 160f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 161f9fd7fdbSMatthew G. Knepley { 162f9fd7fdbSMatthew G. Knepley PetscInt q, d; 163f9fd7fdbSMatthew G. Knepley PetscErrorCode ierr; 164f9fd7fdbSMatthew G. Knepley 165f9fd7fdbSMatthew G. Knepley PetscFunctionBegin; 16698c3331eSBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr); 16721454ff5SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n (", quad->numPoints);CHKERRQ(ierr); 16821454ff5SMatthew G. Knepley for (q = 0; q < quad->numPoints; ++q) { 16921454ff5SMatthew G. Knepley for (d = 0; d < quad->dim; ++d) { 170f9fd7fdbSMatthew G. Knepley if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 171ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 172f9fd7fdbSMatthew G. Knepley } 173ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr); 174f9fd7fdbSMatthew G. Knepley } 175bfa639d9SMatthew G. Knepley PetscFunctionReturn(0); 176bfa639d9SMatthew G. Knepley } 177bfa639d9SMatthew G. Knepley 178bfa639d9SMatthew G. Knepley #undef __FUNCT__ 17937045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval" 18037045ce4SJed Brown /*@ 18137045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 18237045ce4SJed Brown 18337045ce4SJed Brown Not Collective 18437045ce4SJed Brown 18537045ce4SJed Brown Input Arguments: 18637045ce4SJed Brown + npoints - number of spatial points to evaluate at 18737045ce4SJed Brown . points - array of locations to evaluate at 18837045ce4SJed Brown . ndegree - number of basis degrees to evaluate 18937045ce4SJed Brown - degrees - sorted array of degrees to evaluate 19037045ce4SJed Brown 19137045ce4SJed Brown Output Arguments: 1920298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 1930298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 1940298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 19537045ce4SJed Brown 19637045ce4SJed Brown Level: intermediate 19737045ce4SJed Brown 19837045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 19937045ce4SJed Brown @*/ 20037045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 20137045ce4SJed Brown { 20237045ce4SJed Brown PetscInt i,maxdegree; 20337045ce4SJed Brown 20437045ce4SJed Brown PetscFunctionBegin; 20537045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 20637045ce4SJed Brown maxdegree = degrees[ndegree-1]; 20737045ce4SJed Brown for (i=0; i<npoints; i++) { 20837045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 20937045ce4SJed Brown PetscInt j,k; 21037045ce4SJed Brown x = points[i]; 21137045ce4SJed Brown pm2 = 0; 21237045ce4SJed Brown pm1 = 1; 21337045ce4SJed Brown pd2 = 0; 21437045ce4SJed Brown pd1 = 0; 21537045ce4SJed Brown pdd2 = 0; 21637045ce4SJed Brown pdd1 = 0; 21737045ce4SJed Brown k = 0; 21837045ce4SJed Brown if (degrees[k] == 0) { 21937045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 22037045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 22137045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 22237045ce4SJed Brown k++; 22337045ce4SJed Brown } 22437045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 22537045ce4SJed Brown PetscReal p,d,dd; 22637045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 22737045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 22837045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 22937045ce4SJed Brown pm2 = pm1; 23037045ce4SJed Brown pm1 = p; 23137045ce4SJed Brown pd2 = pd1; 23237045ce4SJed Brown pd1 = d; 23337045ce4SJed Brown pdd2 = pdd1; 23437045ce4SJed Brown pdd1 = dd; 23537045ce4SJed Brown if (degrees[k] == j) { 23637045ce4SJed Brown if (B) B[i*ndegree+k] = p; 23737045ce4SJed Brown if (D) D[i*ndegree+k] = d; 23837045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 23937045ce4SJed Brown } 24037045ce4SJed Brown } 24137045ce4SJed Brown } 24237045ce4SJed Brown PetscFunctionReturn(0); 24337045ce4SJed Brown } 24437045ce4SJed Brown 24537045ce4SJed Brown #undef __FUNCT__ 24637045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature" 24737045ce4SJed Brown /*@ 24837045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 24937045ce4SJed Brown 25037045ce4SJed Brown Not Collective 25137045ce4SJed Brown 25237045ce4SJed Brown Input Arguments: 25337045ce4SJed Brown + npoints - number of points 25437045ce4SJed Brown . a - left end of interval (often-1) 25537045ce4SJed Brown - b - right end of interval (often +1) 25637045ce4SJed Brown 25737045ce4SJed Brown Output Arguments: 25837045ce4SJed Brown + x - quadrature points 25937045ce4SJed Brown - w - quadrature weights 26037045ce4SJed Brown 26137045ce4SJed Brown Level: intermediate 26237045ce4SJed Brown 26337045ce4SJed Brown References: 26437045ce4SJed Brown Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 26537045ce4SJed Brown 26637045ce4SJed Brown .seealso: PetscDTLegendreEval() 26737045ce4SJed Brown @*/ 26837045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 26937045ce4SJed Brown { 27037045ce4SJed Brown PetscErrorCode ierr; 27137045ce4SJed Brown PetscInt i; 27237045ce4SJed Brown PetscReal *work; 27337045ce4SJed Brown PetscScalar *Z; 27437045ce4SJed Brown PetscBLASInt N,LDZ,info; 27537045ce4SJed Brown 27637045ce4SJed Brown PetscFunctionBegin; 2770bfcf5a5SMatthew G. Knepley ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 27837045ce4SJed Brown /* Set up the Golub-Welsch system */ 27937045ce4SJed Brown for (i=0; i<npoints; i++) { 28037045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 28137045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 28237045ce4SJed Brown } 283dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 284c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 28537045ce4SJed Brown LDZ = N; 28637045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 2878b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 28837045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 2891c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 29037045ce4SJed Brown 29137045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 29237045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 29337045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 29437045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 2950d644c17SKarl Rupp 29688393a60SJed Brown w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 29737045ce4SJed Brown } 29837045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 29937045ce4SJed Brown PetscFunctionReturn(0); 30037045ce4SJed Brown } 301194825f6SJed Brown 302194825f6SJed Brown #undef __FUNCT__ 303744bafbcSMatthew G. Knepley #define __FUNCT__ "PetscDTGaussTensorQuadrature" 304744bafbcSMatthew G. Knepley /*@ 305744bafbcSMatthew G. Knepley PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 306744bafbcSMatthew G. Knepley 307744bafbcSMatthew G. Knepley Not Collective 308744bafbcSMatthew G. Knepley 309744bafbcSMatthew G. Knepley Input Arguments: 310744bafbcSMatthew G. Knepley + dim - The spatial dimension 311744bafbcSMatthew G. Knepley . npoints - number of points in one dimension 312744bafbcSMatthew G. Knepley . a - left end of interval (often-1) 313744bafbcSMatthew G. Knepley - b - right end of interval (often +1) 314744bafbcSMatthew G. Knepley 315744bafbcSMatthew G. Knepley Output Argument: 316744bafbcSMatthew G. Knepley . q - A PetscQuadrature object 317744bafbcSMatthew G. Knepley 318744bafbcSMatthew G. Knepley Level: intermediate 319744bafbcSMatthew G. Knepley 320744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 321744bafbcSMatthew G. Knepley @*/ 322744bafbcSMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 323744bafbcSMatthew G. Knepley { 324744bafbcSMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k; 325744bafbcSMatthew G. Knepley PetscReal *x, *w, *xw, *ww; 326744bafbcSMatthew G. Knepley PetscErrorCode ierr; 327744bafbcSMatthew G. Knepley 328744bafbcSMatthew G. Knepley PetscFunctionBegin; 329744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 330744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr); 331744bafbcSMatthew G. Knepley /* Set up the Golub-Welsch system */ 332744bafbcSMatthew G. Knepley switch (dim) { 333744bafbcSMatthew G. Knepley case 0: 334744bafbcSMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 335744bafbcSMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 336744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 337744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 338744bafbcSMatthew G. Knepley x[0] = 0.0; 339744bafbcSMatthew G. Knepley w[0] = 1.0; 340744bafbcSMatthew G. Knepley break; 341744bafbcSMatthew G. Knepley case 1: 342744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr); 343744bafbcSMatthew G. Knepley break; 344744bafbcSMatthew G. Knepley case 2: 345744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 346744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 347744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 348744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 349744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+0] = xw[i]; 350744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+1] = xw[j]; 351744bafbcSMatthew G. Knepley w[i*npoints+j] = ww[i] * ww[j]; 352744bafbcSMatthew G. Knepley } 353744bafbcSMatthew G. Knepley } 354744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 355744bafbcSMatthew G. Knepley break; 356744bafbcSMatthew G. Knepley case 3: 357744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 358744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 359744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 360744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 361744bafbcSMatthew G. Knepley for (k = 0; k < npoints; ++k) { 362744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 363744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 364744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 365744bafbcSMatthew G. Knepley w[(i*npoints+j)*npoints+k] = ww[i] * ww[j] * ww[k]; 366744bafbcSMatthew G. Knepley } 367744bafbcSMatthew G. Knepley } 368744bafbcSMatthew G. Knepley } 369744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 370744bafbcSMatthew G. Knepley break; 371744bafbcSMatthew G. Knepley default: 372744bafbcSMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 373744bafbcSMatthew G. Knepley } 374744bafbcSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 375*bcede257SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*q, npoints);CHKERRQ(ierr); 376744bafbcSMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr); 377744bafbcSMatthew G. Knepley PetscFunctionReturn(0); 378744bafbcSMatthew G. Knepley } 379744bafbcSMatthew G. Knepley 380744bafbcSMatthew G. Knepley #undef __FUNCT__ 381494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal" 382494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 383494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 384494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 385494e7359SMatthew G. Knepley { 386494e7359SMatthew G. Knepley PetscReal f = 1.0; 387494e7359SMatthew G. Knepley PetscInt i; 388494e7359SMatthew G. Knepley 389494e7359SMatthew G. Knepley PetscFunctionBegin; 390494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 391494e7359SMatthew G. Knepley *factorial = f; 392494e7359SMatthew G. Knepley PetscFunctionReturn(0); 393494e7359SMatthew G. Knepley } 394494e7359SMatthew G. Knepley 395494e7359SMatthew G. Knepley #undef __FUNCT__ 396494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi" 397494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 398494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 399494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 400494e7359SMatthew G. Knepley { 401494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 402494e7359SMatthew G. Knepley PetscInt k; 403494e7359SMatthew G. Knepley 404494e7359SMatthew G. Knepley PetscFunctionBegin; 405494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 406494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 407494e7359SMatthew G. Knepley apb = a + b; 408494e7359SMatthew G. Knepley pn2 = 1.0; 409494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 410494e7359SMatthew G. Knepley *P = 0.0; 411494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 412494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 413494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 414494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 415494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 416494e7359SMatthew G. Knepley 417494e7359SMatthew G. Knepley a2 = a2 / a1; 418494e7359SMatthew G. Knepley a3 = a3 / a1; 419494e7359SMatthew G. Knepley a4 = a4 / a1; 420494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 421494e7359SMatthew G. Knepley pn2 = pn1; 422494e7359SMatthew G. Knepley pn1 = *P; 423494e7359SMatthew G. Knepley } 424494e7359SMatthew G. Knepley PetscFunctionReturn(0); 425494e7359SMatthew G. Knepley } 426494e7359SMatthew G. Knepley 427494e7359SMatthew G. Knepley #undef __FUNCT__ 428494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative" 429494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 430494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 431494e7359SMatthew G. Knepley { 432494e7359SMatthew G. Knepley PetscReal nP; 433494e7359SMatthew G. Knepley PetscErrorCode ierr; 434494e7359SMatthew G. Knepley 435494e7359SMatthew G. Knepley PetscFunctionBegin; 436494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 437494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 438494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 439494e7359SMatthew G. Knepley PetscFunctionReturn(0); 440494e7359SMatthew G. Knepley } 441494e7359SMatthew G. Knepley 442494e7359SMatthew G. Knepley #undef __FUNCT__ 443494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 444494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 445494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 446494e7359SMatthew G. Knepley { 447494e7359SMatthew G. Knepley PetscFunctionBegin; 448494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 449494e7359SMatthew G. Knepley *eta = y; 450494e7359SMatthew G. Knepley PetscFunctionReturn(0); 451494e7359SMatthew G. Knepley } 452494e7359SMatthew G. Knepley 453494e7359SMatthew G. Knepley #undef __FUNCT__ 454494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 455494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 456494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 457494e7359SMatthew G. Knepley { 458494e7359SMatthew G. Knepley PetscFunctionBegin; 459494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 460494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 461494e7359SMatthew G. Knepley *zeta = z; 462494e7359SMatthew G. Knepley PetscFunctionReturn(0); 463494e7359SMatthew G. Knepley } 464494e7359SMatthew G. Knepley 465494e7359SMatthew G. Knepley #undef __FUNCT__ 466494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 467494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 468494e7359SMatthew G. Knepley { 469494e7359SMatthew G. Knepley PetscInt maxIter = 100; 470494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 471a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 472494e7359SMatthew G. Knepley PetscInt k; 473494e7359SMatthew G. Knepley PetscErrorCode ierr; 474494e7359SMatthew G. Knepley 475494e7359SMatthew G. Knepley PetscFunctionBegin; 476a8291ba1SSatish Balay 4778b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 478a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 4790646a658SBarry Smith a2 = PetscTGamma(a + npoints + 1); 4800646a658SBarry Smith a3 = PetscTGamma(b + npoints + 1); 4810646a658SBarry Smith a4 = PetscTGamma(a + b + npoints + 1); 482a8291ba1SSatish Balay #else 483a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 484a8291ba1SSatish Balay #endif 485a8291ba1SSatish Balay 486494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 487494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 488494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 489494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 490494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 4918b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 492494e7359SMatthew G. Knepley PetscInt j; 493494e7359SMatthew G. Knepley 494494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 495494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 496494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 497494e7359SMatthew G. Knepley PetscInt i; 498494e7359SMatthew G. Knepley 499494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 500494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 501494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 502494e7359SMatthew G. Knepley delta = f / (fp - f * s); 503494e7359SMatthew G. Knepley r = r - delta; 50477b4d14cSPeter Brune if (PetscAbsReal(delta) < eps) break; 505494e7359SMatthew G. Knepley } 506494e7359SMatthew G. Knepley x[k] = r; 507494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 508494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 509494e7359SMatthew G. Knepley } 510494e7359SMatthew G. Knepley PetscFunctionReturn(0); 511494e7359SMatthew G. Knepley } 512494e7359SMatthew G. Knepley 513494e7359SMatthew G. Knepley #undef __FUNCT__ 514494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 515fd9d31fbSMatthew G. Knepley /*@C 516494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 517494e7359SMatthew G. Knepley 518494e7359SMatthew G. Knepley Not Collective 519494e7359SMatthew G. Knepley 520494e7359SMatthew G. Knepley Input Arguments: 521494e7359SMatthew G. Knepley + dim - The simplex dimension 522744bafbcSMatthew G. Knepley . order - The number of points in one dimension 523494e7359SMatthew G. Knepley . a - left end of interval (often-1) 524494e7359SMatthew G. Knepley - b - right end of interval (often +1) 525494e7359SMatthew G. Knepley 526744bafbcSMatthew G. Knepley Output Argument: 527552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 528494e7359SMatthew G. Knepley 529494e7359SMatthew G. Knepley Level: intermediate 530494e7359SMatthew G. Knepley 531494e7359SMatthew G. Knepley References: 532494e7359SMatthew G. Knepley Karniadakis and Sherwin. 533494e7359SMatthew G. Knepley FIAT 534494e7359SMatthew G. Knepley 535744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 536494e7359SMatthew G. Knepley @*/ 537552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q) 538494e7359SMatthew G. Knepley { 539552aa4f7SMatthew G. Knepley PetscInt npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order; 540494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 541494e7359SMatthew G. Knepley PetscInt i, j, k; 542494e7359SMatthew G. Knepley PetscErrorCode ierr; 543494e7359SMatthew G. Knepley 544494e7359SMatthew G. Knepley PetscFunctionBegin; 545494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 546785e854fSJed Brown ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 547785e854fSJed Brown ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 548494e7359SMatthew G. Knepley switch (dim) { 549707aa5c5SMatthew G. Knepley case 0: 550707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 551707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 552785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 553785e854fSJed Brown ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 554707aa5c5SMatthew G. Knepley x[0] = 0.0; 555707aa5c5SMatthew G. Knepley w[0] = 1.0; 556707aa5c5SMatthew G. Knepley break; 557494e7359SMatthew G. Knepley case 1: 558552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr); 559494e7359SMatthew G. Knepley break; 560494e7359SMatthew G. Knepley case 2: 561dcca6d9dSJed Brown ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr); 562552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 563552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 564552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 565552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 566552aa4f7SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr); 567552aa4f7SMatthew G. Knepley w[i*order+j] = 0.5 * wx[i] * wy[j]; 568494e7359SMatthew G. Knepley } 569494e7359SMatthew G. Knepley } 570494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 571494e7359SMatthew G. Knepley break; 572494e7359SMatthew G. Knepley case 3: 573dcca6d9dSJed Brown ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr); 574552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 575552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 576552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 577552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 578552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 579552aa4f7SMatthew G. Knepley for (k = 0; k < order; ++k) { 580552aa4f7SMatthew 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); 581552aa4f7SMatthew G. Knepley w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k]; 582494e7359SMatthew G. Knepley } 583494e7359SMatthew G. Knepley } 584494e7359SMatthew G. Knepley } 585494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 586494e7359SMatthew G. Knepley break; 587494e7359SMatthew G. Knepley default: 588494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 589494e7359SMatthew G. Knepley } 59021454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 591*bcede257SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*q, order);CHKERRQ(ierr); 59221454ff5SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr); 593494e7359SMatthew G. Knepley PetscFunctionReturn(0); 594494e7359SMatthew G. Knepley } 595494e7359SMatthew G. Knepley 596494e7359SMatthew G. Knepley #undef __FUNCT__ 597194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR" 598194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 599194825f6SJed Brown * A in column-major format 600194825f6SJed Brown * Ainv in row-major format 601194825f6SJed Brown * tau has length m 602194825f6SJed Brown * worksize must be >= max(1,n) 603194825f6SJed Brown */ 604194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 605194825f6SJed Brown { 606194825f6SJed Brown PetscErrorCode ierr; 607194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 608194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 609194825f6SJed Brown 610194825f6SJed Brown PetscFunctionBegin; 611194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 612194825f6SJed Brown { 613194825f6SJed Brown PetscInt i,j; 614dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 615194825f6SJed Brown for (j=0; j<n; j++) { 616194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 617194825f6SJed Brown } 618194825f6SJed Brown mstride = m; 619194825f6SJed Brown } 620194825f6SJed Brown #else 621194825f6SJed Brown A = A_in; 622194825f6SJed Brown Ainv = Ainv_out; 623194825f6SJed Brown #endif 624194825f6SJed Brown 625194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 626194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 627194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 628194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 629194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 630001a771dSBarry Smith PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 631194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 632194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 633194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 634194825f6SJed Brown 635194825f6SJed Brown /* Extract an explicit representation of Q */ 636194825f6SJed Brown Q = Ainv; 637194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 638194825f6SJed Brown K = N; /* full rank */ 639001a771dSBarry Smith PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 640194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 641194825f6SJed Brown 642194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 643194825f6SJed Brown Alpha = 1.0; 644194825f6SJed Brown ldb = lda; 645001a771dSBarry Smith PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 646194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 647194825f6SJed Brown 648194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 649194825f6SJed Brown { 650194825f6SJed Brown PetscInt i; 651194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 652194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 653194825f6SJed Brown } 654194825f6SJed Brown #endif 655194825f6SJed Brown PetscFunctionReturn(0); 656194825f6SJed Brown } 657194825f6SJed Brown 658194825f6SJed Brown #undef __FUNCT__ 659194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate" 660194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 661194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 662194825f6SJed Brown { 663194825f6SJed Brown PetscErrorCode ierr; 664194825f6SJed Brown PetscReal *Bv; 665194825f6SJed Brown PetscInt i,j; 666194825f6SJed Brown 667194825f6SJed Brown PetscFunctionBegin; 668785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 669194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 670194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 671194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 672194825f6SJed Brown for (i=0; i<ninterval; i++) { 673194825f6SJed Brown for (j=0; j<ndegree; j++) { 674194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 675194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 676194825f6SJed Brown } 677194825f6SJed Brown } 678194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 679194825f6SJed Brown PetscFunctionReturn(0); 680194825f6SJed Brown } 681194825f6SJed Brown 682194825f6SJed Brown #undef __FUNCT__ 683194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly" 684194825f6SJed Brown /*@ 685194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 686194825f6SJed Brown 687194825f6SJed Brown Not Collective 688194825f6SJed Brown 689194825f6SJed Brown Input Arguments: 690194825f6SJed Brown + degree - degree of reconstruction polynomial 691194825f6SJed Brown . nsource - number of source intervals 692194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 693194825f6SJed Brown . ntarget - number of target intervals 694194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 695194825f6SJed Brown 696194825f6SJed Brown Output Arguments: 697194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 698194825f6SJed Brown 699194825f6SJed Brown Level: advanced 700194825f6SJed Brown 701194825f6SJed Brown .seealso: PetscDTLegendreEval() 702194825f6SJed Brown @*/ 703194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 704194825f6SJed Brown { 705194825f6SJed Brown PetscErrorCode ierr; 706194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 707194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 708194825f6SJed Brown PetscScalar *tau,*work; 709194825f6SJed Brown 710194825f6SJed Brown PetscFunctionBegin; 711194825f6SJed Brown PetscValidRealPointer(sourcex,3); 712194825f6SJed Brown PetscValidRealPointer(targetx,5); 713194825f6SJed Brown PetscValidRealPointer(R,6); 714194825f6SJed 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); 715194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 716194825f6SJed Brown for (i=0; i<nsource; i++) { 71757622a8eSBarry 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]); 718194825f6SJed Brown } 719194825f6SJed Brown for (i=0; i<ntarget; i++) { 72057622a8eSBarry 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]); 721194825f6SJed Brown } 722194825f6SJed Brown #endif 723194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 724194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 725194825f6SJed Brown center = (xmin + xmax)/2; 726194825f6SJed Brown hscale = (xmax - xmin)/2; 727194825f6SJed Brown worksize = nsource; 728dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 729dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 730194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 731194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 732194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 733194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 734194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 735194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 736194825f6SJed Brown for (i=0; i<ntarget; i++) { 737194825f6SJed Brown PetscReal rowsum = 0; 738194825f6SJed Brown for (j=0; j<nsource; j++) { 739194825f6SJed Brown PetscReal sum = 0; 740194825f6SJed Brown for (k=0; k<degree+1; k++) { 741194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 742194825f6SJed Brown } 743194825f6SJed Brown R[i*nsource+j] = sum; 744194825f6SJed Brown rowsum += sum; 745194825f6SJed Brown } 746194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 747194825f6SJed Brown } 748194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 749194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 750194825f6SJed Brown PetscFunctionReturn(0); 751194825f6SJed Brown } 752