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" 2840d8ff71SMatthew G. Knepley /*@ 2940d8ff71SMatthew G. Knepley PetscQuadratureCreate - Create a PetscQuadrature object 3040d8ff71SMatthew G. Knepley 3140d8ff71SMatthew G. Knepley Collective on MPI_Comm 3240d8ff71SMatthew G. Knepley 3340d8ff71SMatthew G. Knepley Input Parameter: 3440d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object 3540d8ff71SMatthew G. Knepley 3640d8ff71SMatthew G. Knepley Output Parameter: 3740d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 3840d8ff71SMatthew G. Knepley 3940d8ff71SMatthew G. Knepley Level: beginner 4040d8ff71SMatthew G. Knepley 4140d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, create 4240d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 4340d8ff71SMatthew G. Knepley @*/ 4421454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 4521454ff5SMatthew G. Knepley { 4621454ff5SMatthew G. Knepley PetscErrorCode ierr; 4721454ff5SMatthew G. Knepley 4821454ff5SMatthew G. Knepley PetscFunctionBegin; 4921454ff5SMatthew G. Knepley PetscValidPointer(q, 2); 5021454ff5SMatthew G. Knepley ierr = DMInitializePackage();CHKERRQ(ierr); 5121454ff5SMatthew G. Knepley ierr = PetscHeaderCreate(*q,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 5221454ff5SMatthew G. Knepley (*q)->dim = -1; 53bcede257SMatthew G. Knepley (*q)->order = -1; 5421454ff5SMatthew G. Knepley (*q)->numPoints = 0; 5521454ff5SMatthew G. Knepley (*q)->points = NULL; 5621454ff5SMatthew G. Knepley (*q)->weights = NULL; 5721454ff5SMatthew G. Knepley PetscFunctionReturn(0); 5821454ff5SMatthew G. Knepley } 5921454ff5SMatthew G. Knepley 6021454ff5SMatthew G. Knepley #undef __FUNCT__ 61*c9638911SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDuplicate" 62*c9638911SMatthew G. Knepley /*@ 63*c9638911SMatthew G. Knepley PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 64*c9638911SMatthew G. Knepley 65*c9638911SMatthew G. Knepley Collective on PetscQuadrature 66*c9638911SMatthew G. Knepley 67*c9638911SMatthew G. Knepley Input Parameter: 68*c9638911SMatthew G. Knepley . q - The PetscQuadrature object 69*c9638911SMatthew G. Knepley 70*c9638911SMatthew G. Knepley Output Parameter: 71*c9638911SMatthew G. Knepley . r - The new PetscQuadrature object 72*c9638911SMatthew G. Knepley 73*c9638911SMatthew G. Knepley Level: beginner 74*c9638911SMatthew G. Knepley 75*c9638911SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, clone 76*c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 77*c9638911SMatthew G. Knepley @*/ 78*c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 79*c9638911SMatthew G. Knepley { 80*c9638911SMatthew G. Knepley PetscInt order, dim, Nq; 81*c9638911SMatthew G. Knepley const PetscReal *points, *weights; 82*c9638911SMatthew G. Knepley PetscReal *p, *w; 83*c9638911SMatthew G. Knepley PetscErrorCode ierr; 84*c9638911SMatthew G. Knepley 85*c9638911SMatthew G. Knepley PetscFunctionBegin; 86*c9638911SMatthew G. Knepley PetscValidPointer(q, 2); 87*c9638911SMatthew G. Knepley ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 88*c9638911SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 89*c9638911SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 90*c9638911SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nq, &points, &weights);CHKERRQ(ierr); 91*c9638911SMatthew G. Knepley ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 92*c9638911SMatthew G. Knepley ierr = PetscMalloc1(Nq, &w);CHKERRQ(ierr); 93*c9638911SMatthew G. Knepley ierr = PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));CHKERRQ(ierr); 94*c9638911SMatthew G. Knepley ierr = PetscMemcpy(w, weights, Nq * sizeof(PetscReal));CHKERRQ(ierr); 95*c9638911SMatthew G. Knepley ierr = PetscQuadratureSetData(*r, dim, Nq, p, w);CHKERRQ(ierr); 96*c9638911SMatthew G. Knepley PetscFunctionReturn(0); 97*c9638911SMatthew G. Knepley } 98*c9638911SMatthew G. Knepley 99*c9638911SMatthew G. Knepley #undef __FUNCT__ 100bfa639d9SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDestroy" 10140d8ff71SMatthew G. Knepley /*@ 10240d8ff71SMatthew G. Knepley PetscQuadratureDestroy - Destroys a PetscQuadrature object 10340d8ff71SMatthew G. Knepley 10440d8ff71SMatthew G. Knepley Collective on PetscQuadrature 10540d8ff71SMatthew G. Knepley 10640d8ff71SMatthew G. Knepley Input Parameter: 10740d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 10840d8ff71SMatthew G. Knepley 10940d8ff71SMatthew G. Knepley Level: beginner 11040d8ff71SMatthew G. Knepley 11140d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, destroy 11240d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 11340d8ff71SMatthew G. Knepley @*/ 114bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 115bfa639d9SMatthew G. Knepley { 116bfa639d9SMatthew G. Knepley PetscErrorCode ierr; 117bfa639d9SMatthew G. Knepley 118bfa639d9SMatthew G. Knepley PetscFunctionBegin; 11921454ff5SMatthew G. Knepley if (!*q) PetscFunctionReturn(0); 12021454ff5SMatthew G. Knepley PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 12121454ff5SMatthew G. Knepley if (--((PetscObject)(*q))->refct > 0) { 12221454ff5SMatthew G. Knepley *q = NULL; 12321454ff5SMatthew G. Knepley PetscFunctionReturn(0); 12421454ff5SMatthew G. Knepley } 12521454ff5SMatthew G. Knepley ierr = PetscFree((*q)->points);CHKERRQ(ierr); 12621454ff5SMatthew G. Knepley ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 12721454ff5SMatthew G. Knepley ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 12821454ff5SMatthew G. Knepley PetscFunctionReturn(0); 12921454ff5SMatthew G. Knepley } 13021454ff5SMatthew G. Knepley 13121454ff5SMatthew G. Knepley #undef __FUNCT__ 132bcede257SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetOrder" 133bcede257SMatthew G. Knepley /*@ 134bcede257SMatthew G. Knepley PetscQuadratureGetOrder - Return the quadrature information 135bcede257SMatthew G. Knepley 136bcede257SMatthew G. Knepley Not collective 137bcede257SMatthew G. Knepley 138bcede257SMatthew G. Knepley Input Parameter: 139bcede257SMatthew G. Knepley . q - The PetscQuadrature object 140bcede257SMatthew G. Knepley 141bcede257SMatthew G. Knepley Output Parameter: 142bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 143bcede257SMatthew G. Knepley 144bcede257SMatthew G. Knepley Output Parameter: 145bcede257SMatthew G. Knepley 146bcede257SMatthew G. Knepley Level: intermediate 147bcede257SMatthew G. Knepley 148bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 149bcede257SMatthew G. Knepley @*/ 150bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 151bcede257SMatthew G. Knepley { 152bcede257SMatthew G. Knepley PetscFunctionBegin; 153bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 154bcede257SMatthew G. Knepley PetscValidPointer(order, 2); 155bcede257SMatthew G. Knepley *order = q->order; 156bcede257SMatthew G. Knepley PetscFunctionReturn(0); 157bcede257SMatthew G. Knepley } 158bcede257SMatthew G. Knepley 159bcede257SMatthew G. Knepley #undef __FUNCT__ 160bcede257SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetOrder" 161bcede257SMatthew G. Knepley /*@ 162bcede257SMatthew G. Knepley PetscQuadratureSetOrder - Return the quadrature information 163bcede257SMatthew G. Knepley 164bcede257SMatthew G. Knepley Not collective 165bcede257SMatthew G. Knepley 166bcede257SMatthew G. Knepley Input Parameters: 167bcede257SMatthew G. Knepley + q - The PetscQuadrature object 168bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 169bcede257SMatthew G. Knepley 170bcede257SMatthew G. Knepley Level: intermediate 171bcede257SMatthew G. Knepley 172bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 173bcede257SMatthew G. Knepley @*/ 174bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 175bcede257SMatthew G. Knepley { 176bcede257SMatthew G. Knepley PetscFunctionBegin; 177bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 178bcede257SMatthew G. Knepley q->order = order; 179bcede257SMatthew G. Knepley PetscFunctionReturn(0); 180bcede257SMatthew G. Knepley } 181bcede257SMatthew G. Knepley 182bcede257SMatthew G. Knepley #undef __FUNCT__ 18321454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetData" 18440d8ff71SMatthew G. Knepley /*@C 18540d8ff71SMatthew G. Knepley PetscQuadratureGetData - Returns the data defining the quadrature 18640d8ff71SMatthew G. Knepley 18740d8ff71SMatthew G. Knepley Not collective 18840d8ff71SMatthew G. Knepley 18940d8ff71SMatthew G. Knepley Input Parameter: 19040d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 19140d8ff71SMatthew G. Knepley 19240d8ff71SMatthew G. Knepley Output Parameters: 19340d8ff71SMatthew G. Knepley + dim - The spatial dimension 19440d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 19540d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 19640d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 19740d8ff71SMatthew G. Knepley 19840d8ff71SMatthew G. Knepley Level: intermediate 19940d8ff71SMatthew G. Knepley 20040d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature 20140d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 20240d8ff71SMatthew G. Knepley @*/ 20321454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 20421454ff5SMatthew G. Knepley { 20521454ff5SMatthew G. Knepley PetscFunctionBegin; 20621454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 20721454ff5SMatthew G. Knepley if (dim) { 20821454ff5SMatthew G. Knepley PetscValidPointer(dim, 2); 20921454ff5SMatthew G. Knepley *dim = q->dim; 21021454ff5SMatthew G. Knepley } 21121454ff5SMatthew G. Knepley if (npoints) { 21221454ff5SMatthew G. Knepley PetscValidPointer(npoints, 3); 21321454ff5SMatthew G. Knepley *npoints = q->numPoints; 21421454ff5SMatthew G. Knepley } 21521454ff5SMatthew G. Knepley if (points) { 21621454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 21721454ff5SMatthew G. Knepley *points = q->points; 21821454ff5SMatthew G. Knepley } 21921454ff5SMatthew G. Knepley if (weights) { 22021454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 22121454ff5SMatthew G. Knepley *weights = q->weights; 22221454ff5SMatthew G. Knepley } 22321454ff5SMatthew G. Knepley PetscFunctionReturn(0); 22421454ff5SMatthew G. Knepley } 22521454ff5SMatthew G. Knepley 22621454ff5SMatthew G. Knepley #undef __FUNCT__ 22721454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetData" 22840d8ff71SMatthew G. Knepley /*@C 22940d8ff71SMatthew G. Knepley PetscQuadratureSetData - Sets the data defining the quadrature 23040d8ff71SMatthew G. Knepley 23140d8ff71SMatthew G. Knepley Not collective 23240d8ff71SMatthew G. Knepley 23340d8ff71SMatthew G. Knepley Input Parameters: 23440d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 23540d8ff71SMatthew G. Knepley . dim - The spatial dimension 23640d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 23740d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 23840d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 23940d8ff71SMatthew G. Knepley 24040d8ff71SMatthew G. Knepley Level: intermediate 24140d8ff71SMatthew G. Knepley 24240d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature 24340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 24440d8ff71SMatthew G. Knepley @*/ 24521454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 24621454ff5SMatthew G. Knepley { 24721454ff5SMatthew G. Knepley PetscFunctionBegin; 24821454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 24921454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 25021454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 25121454ff5SMatthew G. Knepley if (points) { 25221454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 25321454ff5SMatthew G. Knepley q->points = points; 25421454ff5SMatthew G. Knepley } 25521454ff5SMatthew G. Knepley if (weights) { 25621454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 25721454ff5SMatthew G. Knepley q->weights = weights; 25821454ff5SMatthew G. Knepley } 259f9fd7fdbSMatthew G. Knepley PetscFunctionReturn(0); 260f9fd7fdbSMatthew G. Knepley } 261f9fd7fdbSMatthew G. Knepley 262f9fd7fdbSMatthew G. Knepley #undef __FUNCT__ 263f9fd7fdbSMatthew G. Knepley #define __FUNCT__ "PetscQuadratureView" 26440d8ff71SMatthew G. Knepley /*@C 26540d8ff71SMatthew G. Knepley PetscQuadratureView - Views a PetscQuadrature object 26640d8ff71SMatthew G. Knepley 26740d8ff71SMatthew G. Knepley Collective on PetscQuadrature 26840d8ff71SMatthew G. Knepley 26940d8ff71SMatthew G. Knepley Input Parameters: 27040d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 27140d8ff71SMatthew G. Knepley - viewer - The PetscViewer object 27240d8ff71SMatthew G. Knepley 27340d8ff71SMatthew G. Knepley Level: beginner 27440d8ff71SMatthew G. Knepley 27540d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, view 27640d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 27740d8ff71SMatthew G. Knepley @*/ 278f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 279f9fd7fdbSMatthew G. Knepley { 280f9fd7fdbSMatthew G. Knepley PetscInt q, d; 281f9fd7fdbSMatthew G. Knepley PetscErrorCode ierr; 282f9fd7fdbSMatthew G. Knepley 283f9fd7fdbSMatthew G. Knepley PetscFunctionBegin; 28498c3331eSBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr); 28521454ff5SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n (", quad->numPoints);CHKERRQ(ierr); 28621454ff5SMatthew G. Knepley for (q = 0; q < quad->numPoints; ++q) { 28721454ff5SMatthew G. Knepley for (d = 0; d < quad->dim; ++d) { 288f9fd7fdbSMatthew G. Knepley if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 289ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 290f9fd7fdbSMatthew G. Knepley } 291ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr); 292f9fd7fdbSMatthew G. Knepley } 293bfa639d9SMatthew G. Knepley PetscFunctionReturn(0); 294bfa639d9SMatthew G. Knepley } 295bfa639d9SMatthew G. Knepley 296bfa639d9SMatthew G. Knepley #undef __FUNCT__ 29737045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval" 29837045ce4SJed Brown /*@ 29937045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 30037045ce4SJed Brown 30137045ce4SJed Brown Not Collective 30237045ce4SJed Brown 30337045ce4SJed Brown Input Arguments: 30437045ce4SJed Brown + npoints - number of spatial points to evaluate at 30537045ce4SJed Brown . points - array of locations to evaluate at 30637045ce4SJed Brown . ndegree - number of basis degrees to evaluate 30737045ce4SJed Brown - degrees - sorted array of degrees to evaluate 30837045ce4SJed Brown 30937045ce4SJed Brown Output Arguments: 3100298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 3110298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 3120298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 31337045ce4SJed Brown 31437045ce4SJed Brown Level: intermediate 31537045ce4SJed Brown 31637045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 31737045ce4SJed Brown @*/ 31837045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 31937045ce4SJed Brown { 32037045ce4SJed Brown PetscInt i,maxdegree; 32137045ce4SJed Brown 32237045ce4SJed Brown PetscFunctionBegin; 32337045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 32437045ce4SJed Brown maxdegree = degrees[ndegree-1]; 32537045ce4SJed Brown for (i=0; i<npoints; i++) { 32637045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 32737045ce4SJed Brown PetscInt j,k; 32837045ce4SJed Brown x = points[i]; 32937045ce4SJed Brown pm2 = 0; 33037045ce4SJed Brown pm1 = 1; 33137045ce4SJed Brown pd2 = 0; 33237045ce4SJed Brown pd1 = 0; 33337045ce4SJed Brown pdd2 = 0; 33437045ce4SJed Brown pdd1 = 0; 33537045ce4SJed Brown k = 0; 33637045ce4SJed Brown if (degrees[k] == 0) { 33737045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 33837045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 33937045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 34037045ce4SJed Brown k++; 34137045ce4SJed Brown } 34237045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 34337045ce4SJed Brown PetscReal p,d,dd; 34437045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 34537045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 34637045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 34737045ce4SJed Brown pm2 = pm1; 34837045ce4SJed Brown pm1 = p; 34937045ce4SJed Brown pd2 = pd1; 35037045ce4SJed Brown pd1 = d; 35137045ce4SJed Brown pdd2 = pdd1; 35237045ce4SJed Brown pdd1 = dd; 35337045ce4SJed Brown if (degrees[k] == j) { 35437045ce4SJed Brown if (B) B[i*ndegree+k] = p; 35537045ce4SJed Brown if (D) D[i*ndegree+k] = d; 35637045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 35737045ce4SJed Brown } 35837045ce4SJed Brown } 35937045ce4SJed Brown } 36037045ce4SJed Brown PetscFunctionReturn(0); 36137045ce4SJed Brown } 36237045ce4SJed Brown 36337045ce4SJed Brown #undef __FUNCT__ 36437045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature" 36537045ce4SJed Brown /*@ 36637045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 36737045ce4SJed Brown 36837045ce4SJed Brown Not Collective 36937045ce4SJed Brown 37037045ce4SJed Brown Input Arguments: 37137045ce4SJed Brown + npoints - number of points 37237045ce4SJed Brown . a - left end of interval (often-1) 37337045ce4SJed Brown - b - right end of interval (often +1) 37437045ce4SJed Brown 37537045ce4SJed Brown Output Arguments: 37637045ce4SJed Brown + x - quadrature points 37737045ce4SJed Brown - w - quadrature weights 37837045ce4SJed Brown 37937045ce4SJed Brown Level: intermediate 38037045ce4SJed Brown 38137045ce4SJed Brown References: 38237045ce4SJed Brown Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 38337045ce4SJed Brown 38437045ce4SJed Brown .seealso: PetscDTLegendreEval() 38537045ce4SJed Brown @*/ 38637045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 38737045ce4SJed Brown { 38837045ce4SJed Brown PetscErrorCode ierr; 38937045ce4SJed Brown PetscInt i; 39037045ce4SJed Brown PetscReal *work; 39137045ce4SJed Brown PetscScalar *Z; 39237045ce4SJed Brown PetscBLASInt N,LDZ,info; 39337045ce4SJed Brown 39437045ce4SJed Brown PetscFunctionBegin; 3950bfcf5a5SMatthew G. Knepley ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 39637045ce4SJed Brown /* Set up the Golub-Welsch system */ 39737045ce4SJed Brown for (i=0; i<npoints; i++) { 39837045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 39937045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 40037045ce4SJed Brown } 401dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 402c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 40337045ce4SJed Brown LDZ = N; 40437045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 4058b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 40637045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 4071c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 40837045ce4SJed Brown 40937045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 41037045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 41137045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 41237045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 4130d644c17SKarl Rupp 41488393a60SJed Brown w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 41537045ce4SJed Brown } 41637045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 41737045ce4SJed Brown PetscFunctionReturn(0); 41837045ce4SJed Brown } 419194825f6SJed Brown 420194825f6SJed Brown #undef __FUNCT__ 421744bafbcSMatthew G. Knepley #define __FUNCT__ "PetscDTGaussTensorQuadrature" 422744bafbcSMatthew G. Knepley /*@ 423744bafbcSMatthew G. Knepley PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 424744bafbcSMatthew G. Knepley 425744bafbcSMatthew G. Knepley Not Collective 426744bafbcSMatthew G. Knepley 427744bafbcSMatthew G. Knepley Input Arguments: 428744bafbcSMatthew G. Knepley + dim - The spatial dimension 429744bafbcSMatthew G. Knepley . npoints - number of points in one dimension 430744bafbcSMatthew G. Knepley . a - left end of interval (often-1) 431744bafbcSMatthew G. Knepley - b - right end of interval (often +1) 432744bafbcSMatthew G. Knepley 433744bafbcSMatthew G. Knepley Output Argument: 434744bafbcSMatthew G. Knepley . q - A PetscQuadrature object 435744bafbcSMatthew G. Knepley 436744bafbcSMatthew G. Knepley Level: intermediate 437744bafbcSMatthew G. Knepley 438744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 439744bafbcSMatthew G. Knepley @*/ 440744bafbcSMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 441744bafbcSMatthew G. Knepley { 442744bafbcSMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k; 443744bafbcSMatthew G. Knepley PetscReal *x, *w, *xw, *ww; 444744bafbcSMatthew G. Knepley PetscErrorCode ierr; 445744bafbcSMatthew G. Knepley 446744bafbcSMatthew G. Knepley PetscFunctionBegin; 447744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 448744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr); 449744bafbcSMatthew G. Knepley /* Set up the Golub-Welsch system */ 450744bafbcSMatthew G. Knepley switch (dim) { 451744bafbcSMatthew G. Knepley case 0: 452744bafbcSMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 453744bafbcSMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 454744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 455744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 456744bafbcSMatthew G. Knepley x[0] = 0.0; 457744bafbcSMatthew G. Knepley w[0] = 1.0; 458744bafbcSMatthew G. Knepley break; 459744bafbcSMatthew G. Knepley case 1: 460744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr); 461744bafbcSMatthew G. Knepley break; 462744bafbcSMatthew G. Knepley case 2: 463744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 464744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 465744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 466744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 467744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+0] = xw[i]; 468744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+1] = xw[j]; 469744bafbcSMatthew G. Knepley w[i*npoints+j] = ww[i] * ww[j]; 470744bafbcSMatthew G. Knepley } 471744bafbcSMatthew G. Knepley } 472744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 473744bafbcSMatthew G. Knepley break; 474744bafbcSMatthew G. Knepley case 3: 475744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 476744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 477744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 478744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 479744bafbcSMatthew G. Knepley for (k = 0; k < npoints; ++k) { 480744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 481744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 482744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 483744bafbcSMatthew G. Knepley w[(i*npoints+j)*npoints+k] = ww[i] * ww[j] * ww[k]; 484744bafbcSMatthew G. Knepley } 485744bafbcSMatthew G. Knepley } 486744bafbcSMatthew G. Knepley } 487744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 488744bafbcSMatthew G. Knepley break; 489744bafbcSMatthew G. Knepley default: 490744bafbcSMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 491744bafbcSMatthew G. Knepley } 492744bafbcSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 493bcede257SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*q, npoints);CHKERRQ(ierr); 494744bafbcSMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr); 495744bafbcSMatthew G. Knepley PetscFunctionReturn(0); 496744bafbcSMatthew G. Knepley } 497744bafbcSMatthew G. Knepley 498744bafbcSMatthew G. Knepley #undef __FUNCT__ 499494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal" 500494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 501494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 502494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 503494e7359SMatthew G. Knepley { 504494e7359SMatthew G. Knepley PetscReal f = 1.0; 505494e7359SMatthew G. Knepley PetscInt i; 506494e7359SMatthew G. Knepley 507494e7359SMatthew G. Knepley PetscFunctionBegin; 508494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 509494e7359SMatthew G. Knepley *factorial = f; 510494e7359SMatthew G. Knepley PetscFunctionReturn(0); 511494e7359SMatthew G. Knepley } 512494e7359SMatthew G. Knepley 513494e7359SMatthew G. Knepley #undef __FUNCT__ 514494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi" 515494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 516494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 517494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 518494e7359SMatthew G. Knepley { 519494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 520494e7359SMatthew G. Knepley PetscInt k; 521494e7359SMatthew G. Knepley 522494e7359SMatthew G. Knepley PetscFunctionBegin; 523494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 524494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 525494e7359SMatthew G. Knepley apb = a + b; 526494e7359SMatthew G. Knepley pn2 = 1.0; 527494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 528494e7359SMatthew G. Knepley *P = 0.0; 529494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 530494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 531494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 532494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 533494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 534494e7359SMatthew G. Knepley 535494e7359SMatthew G. Knepley a2 = a2 / a1; 536494e7359SMatthew G. Knepley a3 = a3 / a1; 537494e7359SMatthew G. Knepley a4 = a4 / a1; 538494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 539494e7359SMatthew G. Knepley pn2 = pn1; 540494e7359SMatthew G. Knepley pn1 = *P; 541494e7359SMatthew G. Knepley } 542494e7359SMatthew G. Knepley PetscFunctionReturn(0); 543494e7359SMatthew G. Knepley } 544494e7359SMatthew G. Knepley 545494e7359SMatthew G. Knepley #undef __FUNCT__ 546494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative" 547494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 548494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 549494e7359SMatthew G. Knepley { 550494e7359SMatthew G. Knepley PetscReal nP; 551494e7359SMatthew G. Knepley PetscErrorCode ierr; 552494e7359SMatthew G. Knepley 553494e7359SMatthew G. Knepley PetscFunctionBegin; 554494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 555494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 556494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 557494e7359SMatthew G. Knepley PetscFunctionReturn(0); 558494e7359SMatthew G. Knepley } 559494e7359SMatthew G. Knepley 560494e7359SMatthew G. Knepley #undef __FUNCT__ 561494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 562494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 563494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 564494e7359SMatthew G. Knepley { 565494e7359SMatthew G. Knepley PetscFunctionBegin; 566494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 567494e7359SMatthew G. Knepley *eta = y; 568494e7359SMatthew G. Knepley PetscFunctionReturn(0); 569494e7359SMatthew G. Knepley } 570494e7359SMatthew G. Knepley 571494e7359SMatthew G. Knepley #undef __FUNCT__ 572494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 573494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 574494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 575494e7359SMatthew G. Knepley { 576494e7359SMatthew G. Knepley PetscFunctionBegin; 577494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 578494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 579494e7359SMatthew G. Knepley *zeta = z; 580494e7359SMatthew G. Knepley PetscFunctionReturn(0); 581494e7359SMatthew G. Knepley } 582494e7359SMatthew G. Knepley 583494e7359SMatthew G. Knepley #undef __FUNCT__ 584494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 585494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 586494e7359SMatthew G. Knepley { 587494e7359SMatthew G. Knepley PetscInt maxIter = 100; 588494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 589a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 590494e7359SMatthew G. Knepley PetscInt k; 591494e7359SMatthew G. Knepley PetscErrorCode ierr; 592494e7359SMatthew G. Knepley 593494e7359SMatthew G. Knepley PetscFunctionBegin; 594a8291ba1SSatish Balay 5958b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 596a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 5970646a658SBarry Smith a2 = PetscTGamma(a + npoints + 1); 5980646a658SBarry Smith a3 = PetscTGamma(b + npoints + 1); 5990646a658SBarry Smith a4 = PetscTGamma(a + b + npoints + 1); 600a8291ba1SSatish Balay #else 601a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 602a8291ba1SSatish Balay #endif 603a8291ba1SSatish Balay 604494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 605494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 606494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 607494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 608494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 6098b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 610494e7359SMatthew G. Knepley PetscInt j; 611494e7359SMatthew G. Knepley 612494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 613494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 614494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 615494e7359SMatthew G. Knepley PetscInt i; 616494e7359SMatthew G. Knepley 617494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 618494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 619494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 620494e7359SMatthew G. Knepley delta = f / (fp - f * s); 621494e7359SMatthew G. Knepley r = r - delta; 62277b4d14cSPeter Brune if (PetscAbsReal(delta) < eps) break; 623494e7359SMatthew G. Knepley } 624494e7359SMatthew G. Knepley x[k] = r; 625494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 626494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 627494e7359SMatthew G. Knepley } 628494e7359SMatthew G. Knepley PetscFunctionReturn(0); 629494e7359SMatthew G. Knepley } 630494e7359SMatthew G. Knepley 631494e7359SMatthew G. Knepley #undef __FUNCT__ 632494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 633fd9d31fbSMatthew G. Knepley /*@C 634494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 635494e7359SMatthew G. Knepley 636494e7359SMatthew G. Knepley Not Collective 637494e7359SMatthew G. Knepley 638494e7359SMatthew G. Knepley Input Arguments: 639494e7359SMatthew G. Knepley + dim - The simplex dimension 640744bafbcSMatthew G. Knepley . order - The number of points in one dimension 641494e7359SMatthew G. Knepley . a - left end of interval (often-1) 642494e7359SMatthew G. Knepley - b - right end of interval (often +1) 643494e7359SMatthew G. Knepley 644744bafbcSMatthew G. Knepley Output Argument: 645552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 646494e7359SMatthew G. Knepley 647494e7359SMatthew G. Knepley Level: intermediate 648494e7359SMatthew G. Knepley 649494e7359SMatthew G. Knepley References: 650494e7359SMatthew G. Knepley Karniadakis and Sherwin. 651494e7359SMatthew G. Knepley FIAT 652494e7359SMatthew G. Knepley 653744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 654494e7359SMatthew G. Knepley @*/ 655552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q) 656494e7359SMatthew G. Knepley { 657552aa4f7SMatthew G. Knepley PetscInt npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order; 658494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 659494e7359SMatthew G. Knepley PetscInt i, j, k; 660494e7359SMatthew G. Knepley PetscErrorCode ierr; 661494e7359SMatthew G. Knepley 662494e7359SMatthew G. Knepley PetscFunctionBegin; 663494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 664785e854fSJed Brown ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 665785e854fSJed Brown ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 666494e7359SMatthew G. Knepley switch (dim) { 667707aa5c5SMatthew G. Knepley case 0: 668707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 669707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 670785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 671785e854fSJed Brown ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 672707aa5c5SMatthew G. Knepley x[0] = 0.0; 673707aa5c5SMatthew G. Knepley w[0] = 1.0; 674707aa5c5SMatthew G. Knepley break; 675494e7359SMatthew G. Knepley case 1: 676552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr); 677494e7359SMatthew G. Knepley break; 678494e7359SMatthew G. Knepley case 2: 679dcca6d9dSJed Brown ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr); 680552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 681552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 682552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 683552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 684552aa4f7SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr); 685552aa4f7SMatthew G. Knepley w[i*order+j] = 0.5 * wx[i] * wy[j]; 686494e7359SMatthew G. Knepley } 687494e7359SMatthew G. Knepley } 688494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 689494e7359SMatthew G. Knepley break; 690494e7359SMatthew G. Knepley case 3: 691dcca6d9dSJed Brown ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr); 692552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 693552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 694552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 695552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 696552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 697552aa4f7SMatthew G. Knepley for (k = 0; k < order; ++k) { 698552aa4f7SMatthew 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); 699552aa4f7SMatthew G. Knepley w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k]; 700494e7359SMatthew G. Knepley } 701494e7359SMatthew G. Knepley } 702494e7359SMatthew G. Knepley } 703494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 704494e7359SMatthew G. Knepley break; 705494e7359SMatthew G. Knepley default: 706494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 707494e7359SMatthew G. Knepley } 70821454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 709bcede257SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*q, order);CHKERRQ(ierr); 71021454ff5SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr); 711494e7359SMatthew G. Knepley PetscFunctionReturn(0); 712494e7359SMatthew G. Knepley } 713494e7359SMatthew G. Knepley 714494e7359SMatthew G. Knepley #undef __FUNCT__ 715194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR" 716194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 717194825f6SJed Brown * A in column-major format 718194825f6SJed Brown * Ainv in row-major format 719194825f6SJed Brown * tau has length m 720194825f6SJed Brown * worksize must be >= max(1,n) 721194825f6SJed Brown */ 722194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 723194825f6SJed Brown { 724194825f6SJed Brown PetscErrorCode ierr; 725194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 726194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 727194825f6SJed Brown 728194825f6SJed Brown PetscFunctionBegin; 729194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 730194825f6SJed Brown { 731194825f6SJed Brown PetscInt i,j; 732dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 733194825f6SJed Brown for (j=0; j<n; j++) { 734194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 735194825f6SJed Brown } 736194825f6SJed Brown mstride = m; 737194825f6SJed Brown } 738194825f6SJed Brown #else 739194825f6SJed Brown A = A_in; 740194825f6SJed Brown Ainv = Ainv_out; 741194825f6SJed Brown #endif 742194825f6SJed Brown 743194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 744194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 745194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 746194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 747194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 748001a771dSBarry Smith PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 749194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 750194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 751194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 752194825f6SJed Brown 753194825f6SJed Brown /* Extract an explicit representation of Q */ 754194825f6SJed Brown Q = Ainv; 755194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 756194825f6SJed Brown K = N; /* full rank */ 757001a771dSBarry Smith PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 758194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 759194825f6SJed Brown 760194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 761194825f6SJed Brown Alpha = 1.0; 762194825f6SJed Brown ldb = lda; 763001a771dSBarry Smith PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 764194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 765194825f6SJed Brown 766194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 767194825f6SJed Brown { 768194825f6SJed Brown PetscInt i; 769194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 770194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 771194825f6SJed Brown } 772194825f6SJed Brown #endif 773194825f6SJed Brown PetscFunctionReturn(0); 774194825f6SJed Brown } 775194825f6SJed Brown 776194825f6SJed Brown #undef __FUNCT__ 777194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate" 778194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 779194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 780194825f6SJed Brown { 781194825f6SJed Brown PetscErrorCode ierr; 782194825f6SJed Brown PetscReal *Bv; 783194825f6SJed Brown PetscInt i,j; 784194825f6SJed Brown 785194825f6SJed Brown PetscFunctionBegin; 786785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 787194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 788194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 789194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 790194825f6SJed Brown for (i=0; i<ninterval; i++) { 791194825f6SJed Brown for (j=0; j<ndegree; j++) { 792194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 793194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 794194825f6SJed Brown } 795194825f6SJed Brown } 796194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 797194825f6SJed Brown PetscFunctionReturn(0); 798194825f6SJed Brown } 799194825f6SJed Brown 800194825f6SJed Brown #undef __FUNCT__ 801194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly" 802194825f6SJed Brown /*@ 803194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 804194825f6SJed Brown 805194825f6SJed Brown Not Collective 806194825f6SJed Brown 807194825f6SJed Brown Input Arguments: 808194825f6SJed Brown + degree - degree of reconstruction polynomial 809194825f6SJed Brown . nsource - number of source intervals 810194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 811194825f6SJed Brown . ntarget - number of target intervals 812194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 813194825f6SJed Brown 814194825f6SJed Brown Output Arguments: 815194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 816194825f6SJed Brown 817194825f6SJed Brown Level: advanced 818194825f6SJed Brown 819194825f6SJed Brown .seealso: PetscDTLegendreEval() 820194825f6SJed Brown @*/ 821194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 822194825f6SJed Brown { 823194825f6SJed Brown PetscErrorCode ierr; 824194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 825194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 826194825f6SJed Brown PetscScalar *tau,*work; 827194825f6SJed Brown 828194825f6SJed Brown PetscFunctionBegin; 829194825f6SJed Brown PetscValidRealPointer(sourcex,3); 830194825f6SJed Brown PetscValidRealPointer(targetx,5); 831194825f6SJed Brown PetscValidRealPointer(R,6); 832194825f6SJed 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); 833194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 834194825f6SJed Brown for (i=0; i<nsource; i++) { 83557622a8eSBarry 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]); 836194825f6SJed Brown } 837194825f6SJed Brown for (i=0; i<ntarget; i++) { 83857622a8eSBarry 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]); 839194825f6SJed Brown } 840194825f6SJed Brown #endif 841194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 842194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 843194825f6SJed Brown center = (xmin + xmax)/2; 844194825f6SJed Brown hscale = (xmax - xmin)/2; 845194825f6SJed Brown worksize = nsource; 846dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 847dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 848194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 849194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 850194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 851194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 852194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 853194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 854194825f6SJed Brown for (i=0; i<ntarget; i++) { 855194825f6SJed Brown PetscReal rowsum = 0; 856194825f6SJed Brown for (j=0; j<nsource; j++) { 857194825f6SJed Brown PetscReal sum = 0; 858194825f6SJed Brown for (k=0; k<degree+1; k++) { 859194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 860194825f6SJed Brown } 861194825f6SJed Brown R[i*nsource+j] = sum; 862194825f6SJed Brown rowsum += sum; 863194825f6SJed Brown } 864194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 865194825f6SJed Brown } 866194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 867194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 868194825f6SJed Brown PetscFunctionReturn(0); 869194825f6SJed Brown } 870