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 729f144ccSMatthew G. Knepley #ifdef PETSC_HAVE_MPFR 829f144ccSMatthew G. Knepley #include <mpfr.h> 929f144ccSMatthew G. Knepley #endif 10a6fc04d9SSatish Balay 110c35b76eSJed Brown #include <petscdt.h> /*I "petscdt.h" I*/ 1237045ce4SJed Brown #include <petscblaslapack.h> 13af0996ceSBarry Smith #include <petsc/private/petscimpl.h> 14af0996ceSBarry Smith #include <petsc/private/dtimpl.h> 15665c2dedSJed Brown #include <petscviewer.h> 1659804f93SMatthew G. Knepley #include <petscdmplex.h> 1759804f93SMatthew G. Knepley #include <petscdmshell.h> 1837045ce4SJed Brown 190bfcf5a5SMatthew G. Knepley static PetscBool GaussCite = PETSC_FALSE; 200bfcf5a5SMatthew G. Knepley const char GaussCitation[] = "@article{GolubWelsch1969,\n" 210bfcf5a5SMatthew G. Knepley " author = {Golub and Welsch},\n" 220bfcf5a5SMatthew G. Knepley " title = {Calculation of Quadrature Rules},\n" 230bfcf5a5SMatthew G. Knepley " journal = {Math. Comp.},\n" 240bfcf5a5SMatthew G. Knepley " volume = {23},\n" 250bfcf5a5SMatthew G. Knepley " number = {106},\n" 260bfcf5a5SMatthew G. Knepley " pages = {221--230},\n" 270bfcf5a5SMatthew G. Knepley " year = {1969}\n}\n"; 280bfcf5a5SMatthew G. Knepley 2937045ce4SJed Brown #undef __FUNCT__ 3021454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureCreate" 3140d8ff71SMatthew G. Knepley /*@ 3240d8ff71SMatthew G. Knepley PetscQuadratureCreate - Create a PetscQuadrature object 3340d8ff71SMatthew G. Knepley 3440d8ff71SMatthew G. Knepley Collective on MPI_Comm 3540d8ff71SMatthew G. Knepley 3640d8ff71SMatthew G. Knepley Input Parameter: 3740d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object 3840d8ff71SMatthew G. Knepley 3940d8ff71SMatthew G. Knepley Output Parameter: 4040d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 4140d8ff71SMatthew G. Knepley 4240d8ff71SMatthew G. Knepley Level: beginner 4340d8ff71SMatthew G. Knepley 4440d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, create 4540d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 4640d8ff71SMatthew G. Knepley @*/ 4721454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 4821454ff5SMatthew G. Knepley { 4921454ff5SMatthew G. Knepley PetscErrorCode ierr; 5021454ff5SMatthew G. Knepley 5121454ff5SMatthew G. Knepley PetscFunctionBegin; 5221454ff5SMatthew G. Knepley PetscValidPointer(q, 2); 53623436dcSMatthew G. Knepley ierr = PetscSysInitializePackage();CHKERRQ(ierr); 5473107ff1SLisandro Dalcin ierr = PetscHeaderCreate(*q,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 5521454ff5SMatthew G. Knepley (*q)->dim = -1; 56bcede257SMatthew G. Knepley (*q)->order = -1; 5721454ff5SMatthew G. Knepley (*q)->numPoints = 0; 5821454ff5SMatthew G. Knepley (*q)->points = NULL; 5921454ff5SMatthew G. Knepley (*q)->weights = NULL; 6021454ff5SMatthew G. Knepley PetscFunctionReturn(0); 6121454ff5SMatthew G. Knepley } 6221454ff5SMatthew G. Knepley 6321454ff5SMatthew G. Knepley #undef __FUNCT__ 64c9638911SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDuplicate" 65c9638911SMatthew G. Knepley /*@ 66c9638911SMatthew G. Knepley PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 67c9638911SMatthew G. Knepley 68c9638911SMatthew G. Knepley Collective on PetscQuadrature 69c9638911SMatthew G. Knepley 70c9638911SMatthew G. Knepley Input Parameter: 71c9638911SMatthew G. Knepley . q - The PetscQuadrature object 72c9638911SMatthew G. Knepley 73c9638911SMatthew G. Knepley Output Parameter: 74c9638911SMatthew G. Knepley . r - The new PetscQuadrature object 75c9638911SMatthew G. Knepley 76c9638911SMatthew G. Knepley Level: beginner 77c9638911SMatthew G. Knepley 78c9638911SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, clone 79c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 80c9638911SMatthew G. Knepley @*/ 81c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 82c9638911SMatthew G. Knepley { 83c9638911SMatthew G. Knepley PetscInt order, dim, Nq; 84c9638911SMatthew G. Knepley const PetscReal *points, *weights; 85c9638911SMatthew G. Knepley PetscReal *p, *w; 86c9638911SMatthew G. Knepley PetscErrorCode ierr; 87c9638911SMatthew G. Knepley 88c9638911SMatthew G. Knepley PetscFunctionBegin; 89c9638911SMatthew G. Knepley PetscValidPointer(q, 2); 90c9638911SMatthew G. Knepley ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 91c9638911SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 92c9638911SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 93c9638911SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nq, &points, &weights);CHKERRQ(ierr); 94c9638911SMatthew G. Knepley ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 95c9638911SMatthew G. Knepley ierr = PetscMalloc1(Nq, &w);CHKERRQ(ierr); 96c9638911SMatthew G. Knepley ierr = PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));CHKERRQ(ierr); 97c9638911SMatthew G. Knepley ierr = PetscMemcpy(w, weights, Nq * sizeof(PetscReal));CHKERRQ(ierr); 98c9638911SMatthew G. Knepley ierr = PetscQuadratureSetData(*r, dim, Nq, p, w);CHKERRQ(ierr); 99c9638911SMatthew G. Knepley PetscFunctionReturn(0); 100c9638911SMatthew G. Knepley } 101c9638911SMatthew G. Knepley 102c9638911SMatthew G. Knepley #undef __FUNCT__ 103bfa639d9SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDestroy" 10440d8ff71SMatthew G. Knepley /*@ 10540d8ff71SMatthew G. Knepley PetscQuadratureDestroy - Destroys a PetscQuadrature object 10640d8ff71SMatthew G. Knepley 10740d8ff71SMatthew G. Knepley Collective on PetscQuadrature 10840d8ff71SMatthew G. Knepley 10940d8ff71SMatthew G. Knepley Input Parameter: 11040d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 11140d8ff71SMatthew G. Knepley 11240d8ff71SMatthew G. Knepley Level: beginner 11340d8ff71SMatthew G. Knepley 11440d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, destroy 11540d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 11640d8ff71SMatthew G. Knepley @*/ 117bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 118bfa639d9SMatthew G. Knepley { 119bfa639d9SMatthew G. Knepley PetscErrorCode ierr; 120bfa639d9SMatthew G. Knepley 121bfa639d9SMatthew G. Knepley PetscFunctionBegin; 12221454ff5SMatthew G. Knepley if (!*q) PetscFunctionReturn(0); 12321454ff5SMatthew G. Knepley PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 12421454ff5SMatthew G. Knepley if (--((PetscObject)(*q))->refct > 0) { 12521454ff5SMatthew G. Knepley *q = NULL; 12621454ff5SMatthew G. Knepley PetscFunctionReturn(0); 12721454ff5SMatthew G. Knepley } 12821454ff5SMatthew G. Knepley ierr = PetscFree((*q)->points);CHKERRQ(ierr); 12921454ff5SMatthew G. Knepley ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 13021454ff5SMatthew G. Knepley ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 13121454ff5SMatthew G. Knepley PetscFunctionReturn(0); 13221454ff5SMatthew G. Knepley } 13321454ff5SMatthew G. Knepley 13421454ff5SMatthew G. Knepley #undef __FUNCT__ 135bcede257SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetOrder" 136bcede257SMatthew G. Knepley /*@ 137bcede257SMatthew G. Knepley PetscQuadratureGetOrder - Return the quadrature information 138bcede257SMatthew G. Knepley 139bcede257SMatthew G. Knepley Not collective 140bcede257SMatthew G. Knepley 141bcede257SMatthew G. Knepley Input Parameter: 142bcede257SMatthew G. Knepley . q - The PetscQuadrature object 143bcede257SMatthew G. Knepley 144bcede257SMatthew G. Knepley Output Parameter: 145bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 146bcede257SMatthew G. Knepley 147bcede257SMatthew G. Knepley Output Parameter: 148bcede257SMatthew G. Knepley 149bcede257SMatthew G. Knepley Level: intermediate 150bcede257SMatthew G. Knepley 151bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 152bcede257SMatthew G. Knepley @*/ 153bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 154bcede257SMatthew G. Knepley { 155bcede257SMatthew G. Knepley PetscFunctionBegin; 156bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 157bcede257SMatthew G. Knepley PetscValidPointer(order, 2); 158bcede257SMatthew G. Knepley *order = q->order; 159bcede257SMatthew G. Knepley PetscFunctionReturn(0); 160bcede257SMatthew G. Knepley } 161bcede257SMatthew G. Knepley 162bcede257SMatthew G. Knepley #undef __FUNCT__ 163bcede257SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetOrder" 164bcede257SMatthew G. Knepley /*@ 165bcede257SMatthew G. Knepley PetscQuadratureSetOrder - Return the quadrature information 166bcede257SMatthew G. Knepley 167bcede257SMatthew G. Knepley Not collective 168bcede257SMatthew G. Knepley 169bcede257SMatthew G. Knepley Input Parameters: 170bcede257SMatthew G. Knepley + q - The PetscQuadrature object 171bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 172bcede257SMatthew G. Knepley 173bcede257SMatthew G. Knepley Level: intermediate 174bcede257SMatthew G. Knepley 175bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 176bcede257SMatthew G. Knepley @*/ 177bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 178bcede257SMatthew G. Knepley { 179bcede257SMatthew G. Knepley PetscFunctionBegin; 180bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 181bcede257SMatthew G. Knepley q->order = order; 182bcede257SMatthew G. Knepley PetscFunctionReturn(0); 183bcede257SMatthew G. Knepley } 184bcede257SMatthew G. Knepley 185bcede257SMatthew G. Knepley #undef __FUNCT__ 18621454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureGetData" 18740d8ff71SMatthew G. Knepley /*@C 18840d8ff71SMatthew G. Knepley PetscQuadratureGetData - Returns the data defining the quadrature 18940d8ff71SMatthew G. Knepley 19040d8ff71SMatthew G. Knepley Not collective 19140d8ff71SMatthew G. Knepley 19240d8ff71SMatthew G. Knepley Input Parameter: 19340d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 19440d8ff71SMatthew G. Knepley 19540d8ff71SMatthew G. Knepley Output Parameters: 19640d8ff71SMatthew G. Knepley + dim - The spatial dimension 19740d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 19840d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 19940d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 20040d8ff71SMatthew G. Knepley 20140d8ff71SMatthew G. Knepley Level: intermediate 20240d8ff71SMatthew G. Knepley 20340d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature 20440d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 20540d8ff71SMatthew G. Knepley @*/ 20621454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 20721454ff5SMatthew G. Knepley { 20821454ff5SMatthew G. Knepley PetscFunctionBegin; 20921454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 21021454ff5SMatthew G. Knepley if (dim) { 21121454ff5SMatthew G. Knepley PetscValidPointer(dim, 2); 21221454ff5SMatthew G. Knepley *dim = q->dim; 21321454ff5SMatthew G. Knepley } 21421454ff5SMatthew G. Knepley if (npoints) { 21521454ff5SMatthew G. Knepley PetscValidPointer(npoints, 3); 21621454ff5SMatthew G. Knepley *npoints = q->numPoints; 21721454ff5SMatthew G. Knepley } 21821454ff5SMatthew G. Knepley if (points) { 21921454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 22021454ff5SMatthew G. Knepley *points = q->points; 22121454ff5SMatthew G. Knepley } 22221454ff5SMatthew G. Knepley if (weights) { 22321454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 22421454ff5SMatthew G. Knepley *weights = q->weights; 22521454ff5SMatthew G. Knepley } 22621454ff5SMatthew G. Knepley PetscFunctionReturn(0); 22721454ff5SMatthew G. Knepley } 22821454ff5SMatthew G. Knepley 22921454ff5SMatthew G. Knepley #undef __FUNCT__ 23021454ff5SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureSetData" 23140d8ff71SMatthew G. Knepley /*@C 23240d8ff71SMatthew G. Knepley PetscQuadratureSetData - Sets the data defining the quadrature 23340d8ff71SMatthew G. Knepley 23440d8ff71SMatthew G. Knepley Not collective 23540d8ff71SMatthew G. Knepley 23640d8ff71SMatthew G. Knepley Input Parameters: 23740d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 23840d8ff71SMatthew G. Knepley . dim - The spatial dimension 23940d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 24040d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 24140d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 24240d8ff71SMatthew G. Knepley 24340d8ff71SMatthew G. Knepley Level: intermediate 24440d8ff71SMatthew G. Knepley 24540d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature 24640d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 24740d8ff71SMatthew G. Knepley @*/ 24821454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 24921454ff5SMatthew G. Knepley { 25021454ff5SMatthew G. Knepley PetscFunctionBegin; 25121454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 25221454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 25321454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 25421454ff5SMatthew G. Knepley if (points) { 25521454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 25621454ff5SMatthew G. Knepley q->points = points; 25721454ff5SMatthew G. Knepley } 25821454ff5SMatthew G. Knepley if (weights) { 25921454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 26021454ff5SMatthew G. Knepley q->weights = weights; 26121454ff5SMatthew G. Knepley } 262f9fd7fdbSMatthew G. Knepley PetscFunctionReturn(0); 263f9fd7fdbSMatthew G. Knepley } 264f9fd7fdbSMatthew G. Knepley 265f9fd7fdbSMatthew G. Knepley #undef __FUNCT__ 266f9fd7fdbSMatthew G. Knepley #define __FUNCT__ "PetscQuadratureView" 26740d8ff71SMatthew G. Knepley /*@C 26840d8ff71SMatthew G. Knepley PetscQuadratureView - Views a PetscQuadrature object 26940d8ff71SMatthew G. Knepley 27040d8ff71SMatthew G. Knepley Collective on PetscQuadrature 27140d8ff71SMatthew G. Knepley 27240d8ff71SMatthew G. Knepley Input Parameters: 27340d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 27440d8ff71SMatthew G. Knepley - viewer - The PetscViewer object 27540d8ff71SMatthew G. Knepley 27640d8ff71SMatthew G. Knepley Level: beginner 27740d8ff71SMatthew G. Knepley 27840d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, view 27940d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 28040d8ff71SMatthew G. Knepley @*/ 281f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 282f9fd7fdbSMatthew G. Knepley { 283f9fd7fdbSMatthew G. Knepley PetscInt q, d; 284f9fd7fdbSMatthew G. Knepley PetscErrorCode ierr; 285f9fd7fdbSMatthew G. Knepley 286f9fd7fdbSMatthew G. Knepley PetscFunctionBegin; 28798c3331eSBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr); 28821454ff5SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n (", quad->numPoints);CHKERRQ(ierr); 28921454ff5SMatthew G. Knepley for (q = 0; q < quad->numPoints; ++q) { 29021454ff5SMatthew G. Knepley for (d = 0; d < quad->dim; ++d) { 291f9fd7fdbSMatthew G. Knepley if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 292ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 293f9fd7fdbSMatthew G. Knepley } 294ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr); 295f9fd7fdbSMatthew G. Knepley } 296bfa639d9SMatthew G. Knepley PetscFunctionReturn(0); 297bfa639d9SMatthew G. Knepley } 298bfa639d9SMatthew G. Knepley 299bfa639d9SMatthew G. Knepley #undef __FUNCT__ 30089710940SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureExpandComposite" 30189710940SMatthew G. Knepley /*@C 30289710940SMatthew G. Knepley PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 30389710940SMatthew G. Knepley 30489710940SMatthew G. Knepley Not collective 30589710940SMatthew G. Knepley 30689710940SMatthew G. Knepley Input Parameter: 30789710940SMatthew G. Knepley + q - The original PetscQuadrature 30889710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into 30989710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement 31089710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement 31189710940SMatthew G. Knepley 31289710940SMatthew G. Knepley Output Parameters: 31389710940SMatthew G. Knepley . dim - The dimension 31489710940SMatthew G. Knepley 31589710940SMatthew G. Knepley Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 31689710940SMatthew G. Knepley 31789710940SMatthew G. Knepley Level: intermediate 31889710940SMatthew G. Knepley 31989710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 32089710940SMatthew G. Knepley @*/ 32189710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 32289710940SMatthew G. Knepley { 32389710940SMatthew G. Knepley const PetscReal *points, *weights; 32489710940SMatthew G. Knepley PetscReal *pointsRef, *weightsRef; 32589710940SMatthew G. Knepley PetscInt dim, order, npoints, npointsRef, c, p, d, e; 32689710940SMatthew G. Knepley PetscErrorCode ierr; 32789710940SMatthew G. Knepley 32889710940SMatthew G. Knepley PetscFunctionBegin; 32989710940SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 33089710940SMatthew G. Knepley PetscValidPointer(v0, 3); 33189710940SMatthew G. Knepley PetscValidPointer(jac, 4); 33289710940SMatthew G. Knepley PetscValidPointer(qref, 5); 33389710940SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 33489710940SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 33589710940SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &npoints, &points, &weights);CHKERRQ(ierr); 33689710940SMatthew G. Knepley npointsRef = npoints*numSubelements; 33789710940SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 33889710940SMatthew G. Knepley ierr = PetscMalloc1(npointsRef,&weightsRef);CHKERRQ(ierr); 33989710940SMatthew G. Knepley for (c = 0; c < numSubelements; ++c) { 34089710940SMatthew G. Knepley for (p = 0; p < npoints; ++p) { 34189710940SMatthew G. Knepley for (d = 0; d < dim; ++d) { 34289710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 34389710940SMatthew G. Knepley for (e = 0; e < dim; ++e) { 34489710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 34589710940SMatthew G. Knepley } 34689710940SMatthew G. Knepley } 34789710940SMatthew G. Knepley /* Could also use detJ here */ 34889710940SMatthew G. Knepley weightsRef[c*npoints+p] = weights[p]/numSubelements; 34989710940SMatthew G. Knepley } 35089710940SMatthew G. Knepley } 35189710940SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 35289710940SMatthew G. Knepley ierr = PetscQuadratureSetData(*qref, dim, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 35389710940SMatthew G. Knepley PetscFunctionReturn(0); 35489710940SMatthew G. Knepley } 35589710940SMatthew G. Knepley 35689710940SMatthew G. Knepley #undef __FUNCT__ 35737045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval" 35837045ce4SJed Brown /*@ 35937045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 36037045ce4SJed Brown 36137045ce4SJed Brown Not Collective 36237045ce4SJed Brown 36337045ce4SJed Brown Input Arguments: 36437045ce4SJed Brown + npoints - number of spatial points to evaluate at 36537045ce4SJed Brown . points - array of locations to evaluate at 36637045ce4SJed Brown . ndegree - number of basis degrees to evaluate 36737045ce4SJed Brown - degrees - sorted array of degrees to evaluate 36837045ce4SJed Brown 36937045ce4SJed Brown Output Arguments: 3700298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 3710298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 3720298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 37337045ce4SJed Brown 37437045ce4SJed Brown Level: intermediate 37537045ce4SJed Brown 37637045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 37737045ce4SJed Brown @*/ 37837045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 37937045ce4SJed Brown { 38037045ce4SJed Brown PetscInt i,maxdegree; 38137045ce4SJed Brown 38237045ce4SJed Brown PetscFunctionBegin; 38337045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 38437045ce4SJed Brown maxdegree = degrees[ndegree-1]; 38537045ce4SJed Brown for (i=0; i<npoints; i++) { 38637045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 38737045ce4SJed Brown PetscInt j,k; 38837045ce4SJed Brown x = points[i]; 38937045ce4SJed Brown pm2 = 0; 39037045ce4SJed Brown pm1 = 1; 39137045ce4SJed Brown pd2 = 0; 39237045ce4SJed Brown pd1 = 0; 39337045ce4SJed Brown pdd2 = 0; 39437045ce4SJed Brown pdd1 = 0; 39537045ce4SJed Brown k = 0; 39637045ce4SJed Brown if (degrees[k] == 0) { 39737045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 39837045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 39937045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 40037045ce4SJed Brown k++; 40137045ce4SJed Brown } 40237045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 40337045ce4SJed Brown PetscReal p,d,dd; 40437045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 40537045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 40637045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 40737045ce4SJed Brown pm2 = pm1; 40837045ce4SJed Brown pm1 = p; 40937045ce4SJed Brown pd2 = pd1; 41037045ce4SJed Brown pd1 = d; 41137045ce4SJed Brown pdd2 = pdd1; 41237045ce4SJed Brown pdd1 = dd; 41337045ce4SJed Brown if (degrees[k] == j) { 41437045ce4SJed Brown if (B) B[i*ndegree+k] = p; 41537045ce4SJed Brown if (D) D[i*ndegree+k] = d; 41637045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 41737045ce4SJed Brown } 41837045ce4SJed Brown } 41937045ce4SJed Brown } 42037045ce4SJed Brown PetscFunctionReturn(0); 42137045ce4SJed Brown } 42237045ce4SJed Brown 42337045ce4SJed Brown #undef __FUNCT__ 42437045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature" 42537045ce4SJed Brown /*@ 42637045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 42737045ce4SJed Brown 42837045ce4SJed Brown Not Collective 42937045ce4SJed Brown 43037045ce4SJed Brown Input Arguments: 43137045ce4SJed Brown + npoints - number of points 43237045ce4SJed Brown . a - left end of interval (often-1) 43337045ce4SJed Brown - b - right end of interval (often +1) 43437045ce4SJed Brown 43537045ce4SJed Brown Output Arguments: 43637045ce4SJed Brown + x - quadrature points 43737045ce4SJed Brown - w - quadrature weights 43837045ce4SJed Brown 43937045ce4SJed Brown Level: intermediate 44037045ce4SJed Brown 44137045ce4SJed Brown References: 44296a0c994SBarry Smith . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 44337045ce4SJed Brown 44437045ce4SJed Brown .seealso: PetscDTLegendreEval() 44537045ce4SJed Brown @*/ 44637045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 44737045ce4SJed Brown { 44837045ce4SJed Brown PetscErrorCode ierr; 44937045ce4SJed Brown PetscInt i; 45037045ce4SJed Brown PetscReal *work; 45137045ce4SJed Brown PetscScalar *Z; 45237045ce4SJed Brown PetscBLASInt N,LDZ,info; 45337045ce4SJed Brown 45437045ce4SJed Brown PetscFunctionBegin; 4550bfcf5a5SMatthew G. Knepley ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 45637045ce4SJed Brown /* Set up the Golub-Welsch system */ 45737045ce4SJed Brown for (i=0; i<npoints; i++) { 45837045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 45937045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 46037045ce4SJed Brown } 461dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 462c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 46337045ce4SJed Brown LDZ = N; 46437045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 4658b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 46637045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 4671c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 46837045ce4SJed Brown 46937045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 47037045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 47137045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 47219a57d60SBarry Smith if (x[i] == -0.0) x[i] = 0.0; 47337045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 4740d644c17SKarl Rupp 47588393a60SJed Brown w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 47637045ce4SJed Brown } 47737045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 47837045ce4SJed Brown PetscFunctionReturn(0); 47937045ce4SJed Brown } 480194825f6SJed Brown 481194825f6SJed Brown #undef __FUNCT__ 482744bafbcSMatthew G. Knepley #define __FUNCT__ "PetscDTGaussTensorQuadrature" 483744bafbcSMatthew G. Knepley /*@ 484744bafbcSMatthew G. Knepley PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 485744bafbcSMatthew G. Knepley 486744bafbcSMatthew G. Knepley Not Collective 487744bafbcSMatthew G. Knepley 488744bafbcSMatthew G. Knepley Input Arguments: 489744bafbcSMatthew G. Knepley + dim - The spatial dimension 490744bafbcSMatthew G. Knepley . npoints - number of points in one dimension 491744bafbcSMatthew G. Knepley . a - left end of interval (often-1) 492744bafbcSMatthew G. Knepley - b - right end of interval (often +1) 493744bafbcSMatthew G. Knepley 494744bafbcSMatthew G. Knepley Output Argument: 495744bafbcSMatthew G. Knepley . q - A PetscQuadrature object 496744bafbcSMatthew G. Knepley 497744bafbcSMatthew G. Knepley Level: intermediate 498744bafbcSMatthew G. Knepley 499744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 500744bafbcSMatthew G. Knepley @*/ 501744bafbcSMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 502744bafbcSMatthew G. Knepley { 503744bafbcSMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k; 504744bafbcSMatthew G. Knepley PetscReal *x, *w, *xw, *ww; 505744bafbcSMatthew G. Knepley PetscErrorCode ierr; 506744bafbcSMatthew G. Knepley 507744bafbcSMatthew G. Knepley PetscFunctionBegin; 508744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 509744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr); 510744bafbcSMatthew G. Knepley /* Set up the Golub-Welsch system */ 511744bafbcSMatthew G. Knepley switch (dim) { 512744bafbcSMatthew G. Knepley case 0: 513744bafbcSMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 514744bafbcSMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 515744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 516744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 517744bafbcSMatthew G. Knepley x[0] = 0.0; 518744bafbcSMatthew G. Knepley w[0] = 1.0; 519744bafbcSMatthew G. Knepley break; 520744bafbcSMatthew G. Knepley case 1: 521744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr); 522744bafbcSMatthew G. Knepley break; 523744bafbcSMatthew G. Knepley case 2: 524744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 525744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 526744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 527744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 528744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+0] = xw[i]; 529744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+1] = xw[j]; 530744bafbcSMatthew G. Knepley w[i*npoints+j] = ww[i] * ww[j]; 531744bafbcSMatthew G. Knepley } 532744bafbcSMatthew G. Knepley } 533744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 534744bafbcSMatthew G. Knepley break; 535744bafbcSMatthew G. Knepley case 3: 536744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 537744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 538744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 539744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 540744bafbcSMatthew G. Knepley for (k = 0; k < npoints; ++k) { 541744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 542744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 543744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 544744bafbcSMatthew G. Knepley w[(i*npoints+j)*npoints+k] = ww[i] * ww[j] * ww[k]; 545744bafbcSMatthew G. Knepley } 546744bafbcSMatthew G. Knepley } 547744bafbcSMatthew G. Knepley } 548744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 549744bafbcSMatthew G. Knepley break; 550744bafbcSMatthew G. Knepley default: 551744bafbcSMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 552744bafbcSMatthew G. Knepley } 553744bafbcSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 554bcede257SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*q, npoints);CHKERRQ(ierr); 555744bafbcSMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr); 556744bafbcSMatthew G. Knepley PetscFunctionReturn(0); 557744bafbcSMatthew G. Knepley } 558744bafbcSMatthew G. Knepley 559744bafbcSMatthew G. Knepley #undef __FUNCT__ 560494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal" 561494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 562494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 563494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 564494e7359SMatthew G. Knepley { 565494e7359SMatthew G. Knepley PetscReal f = 1.0; 566494e7359SMatthew G. Knepley PetscInt i; 567494e7359SMatthew G. Knepley 568494e7359SMatthew G. Knepley PetscFunctionBegin; 569494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 570494e7359SMatthew G. Knepley *factorial = f; 571494e7359SMatthew G. Knepley PetscFunctionReturn(0); 572494e7359SMatthew G. Knepley } 573494e7359SMatthew G. Knepley 574494e7359SMatthew G. Knepley #undef __FUNCT__ 575494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi" 576494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 577494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 578494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 579494e7359SMatthew G. Knepley { 580494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 581494e7359SMatthew G. Knepley PetscInt k; 582494e7359SMatthew G. Knepley 583494e7359SMatthew G. Knepley PetscFunctionBegin; 584494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 585494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 586494e7359SMatthew G. Knepley apb = a + b; 587494e7359SMatthew G. Knepley pn2 = 1.0; 588494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 589494e7359SMatthew G. Knepley *P = 0.0; 590494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 591494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 592494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 593494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 594494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 595494e7359SMatthew G. Knepley 596494e7359SMatthew G. Knepley a2 = a2 / a1; 597494e7359SMatthew G. Knepley a3 = a3 / a1; 598494e7359SMatthew G. Knepley a4 = a4 / a1; 599494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 600494e7359SMatthew G. Knepley pn2 = pn1; 601494e7359SMatthew G. Knepley pn1 = *P; 602494e7359SMatthew G. Knepley } 603494e7359SMatthew G. Knepley PetscFunctionReturn(0); 604494e7359SMatthew G. Knepley } 605494e7359SMatthew G. Knepley 606494e7359SMatthew G. Knepley #undef __FUNCT__ 607494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative" 608494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 609494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 610494e7359SMatthew G. Knepley { 611494e7359SMatthew G. Knepley PetscReal nP; 612494e7359SMatthew G. Knepley PetscErrorCode ierr; 613494e7359SMatthew G. Knepley 614494e7359SMatthew G. Knepley PetscFunctionBegin; 615494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 616494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 617494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 618494e7359SMatthew G. Knepley PetscFunctionReturn(0); 619494e7359SMatthew G. Knepley } 620494e7359SMatthew G. Knepley 621494e7359SMatthew G. Knepley #undef __FUNCT__ 622494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 623494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 624494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 625494e7359SMatthew G. Knepley { 626494e7359SMatthew G. Knepley PetscFunctionBegin; 627494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 628494e7359SMatthew G. Knepley *eta = y; 629494e7359SMatthew G. Knepley PetscFunctionReturn(0); 630494e7359SMatthew G. Knepley } 631494e7359SMatthew G. Knepley 632494e7359SMatthew G. Knepley #undef __FUNCT__ 633494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 634494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 635494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 636494e7359SMatthew G. Knepley { 637494e7359SMatthew G. Knepley PetscFunctionBegin; 638494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 639494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 640494e7359SMatthew G. Knepley *zeta = z; 641494e7359SMatthew G. Knepley PetscFunctionReturn(0); 642494e7359SMatthew G. Knepley } 643494e7359SMatthew G. Knepley 644494e7359SMatthew G. Knepley #undef __FUNCT__ 645494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 646494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 647494e7359SMatthew G. Knepley { 648494e7359SMatthew G. Knepley PetscInt maxIter = 100; 649494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 650a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 651494e7359SMatthew G. Knepley PetscInt k; 652494e7359SMatthew G. Knepley PetscErrorCode ierr; 653494e7359SMatthew G. Knepley 654494e7359SMatthew G. Knepley PetscFunctionBegin; 655a8291ba1SSatish Balay 6568b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 657a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 6580646a658SBarry Smith a2 = PetscTGamma(a + npoints + 1); 6590646a658SBarry Smith a3 = PetscTGamma(b + npoints + 1); 6600646a658SBarry Smith a4 = PetscTGamma(a + b + npoints + 1); 661a8291ba1SSatish Balay #else 66229bcbfd0SToby Isaac { 663*d24bbb91SToby Isaac PetscInt ia, ib; 66429bcbfd0SToby Isaac 665*d24bbb91SToby Isaac ia = (PetscInt) a; 666*d24bbb91SToby Isaac ib = (PetscInt) b; 667*d24bbb91SToby Isaac if (ia == a && ib == b && ia + npoints + 1 > 0 && ib + npoints + 1 > 0 && ia + ib + npoints + 1 > 0) { /* All gamma(x) terms are (x-1)! terms */ 668*d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr); 669*d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr); 670*d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr); 67129bcbfd0SToby Isaac } else { 672a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 67329bcbfd0SToby Isaac } 67429bcbfd0SToby Isaac } 675a8291ba1SSatish Balay #endif 676a8291ba1SSatish Balay 677494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 678494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 679494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 680494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 681494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 6828b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 683494e7359SMatthew G. Knepley PetscInt j; 684494e7359SMatthew G. Knepley 685494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 686494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 687494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 688494e7359SMatthew G. Knepley PetscInt i; 689494e7359SMatthew G. Knepley 690494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 691494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 692494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 693494e7359SMatthew G. Knepley delta = f / (fp - f * s); 694494e7359SMatthew G. Knepley r = r - delta; 69577b4d14cSPeter Brune if (PetscAbsReal(delta) < eps) break; 696494e7359SMatthew G. Knepley } 697494e7359SMatthew G. Knepley x[k] = r; 698494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 699494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 700494e7359SMatthew G. Knepley } 701494e7359SMatthew G. Knepley PetscFunctionReturn(0); 702494e7359SMatthew G. Knepley } 703494e7359SMatthew G. Knepley 704494e7359SMatthew G. Knepley #undef __FUNCT__ 705494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 706fd9d31fbSMatthew G. Knepley /*@C 707494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 708494e7359SMatthew G. Knepley 709494e7359SMatthew G. Knepley Not Collective 710494e7359SMatthew G. Knepley 711494e7359SMatthew G. Knepley Input Arguments: 712494e7359SMatthew G. Knepley + dim - The simplex dimension 713744bafbcSMatthew G. Knepley . order - The number of points in one dimension 714494e7359SMatthew G. Knepley . a - left end of interval (often-1) 715494e7359SMatthew G. Knepley - b - right end of interval (often +1) 716494e7359SMatthew G. Knepley 717744bafbcSMatthew G. Knepley Output Argument: 718552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 719494e7359SMatthew G. Knepley 720494e7359SMatthew G. Knepley Level: intermediate 721494e7359SMatthew G. Knepley 722494e7359SMatthew G. Knepley References: 72396a0c994SBarry Smith . 1. - Karniadakis and Sherwin. FIAT 724494e7359SMatthew G. Knepley 725744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 726494e7359SMatthew G. Knepley @*/ 727552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q) 728494e7359SMatthew G. Knepley { 729552aa4f7SMatthew G. Knepley PetscInt npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order; 730494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 731494e7359SMatthew G. Knepley PetscInt i, j, k; 732494e7359SMatthew G. Knepley PetscErrorCode ierr; 733494e7359SMatthew G. Knepley 734494e7359SMatthew G. Knepley PetscFunctionBegin; 735494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 736785e854fSJed Brown ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 737785e854fSJed Brown ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 738494e7359SMatthew G. Knepley switch (dim) { 739707aa5c5SMatthew G. Knepley case 0: 740707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 741707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 742785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 743785e854fSJed Brown ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 744707aa5c5SMatthew G. Knepley x[0] = 0.0; 745707aa5c5SMatthew G. Knepley w[0] = 1.0; 746707aa5c5SMatthew G. Knepley break; 747494e7359SMatthew G. Knepley case 1: 748552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr); 749494e7359SMatthew G. Knepley break; 750494e7359SMatthew G. Knepley case 2: 751dcca6d9dSJed Brown ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr); 752552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 753552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 754552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 755552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 756552aa4f7SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr); 757552aa4f7SMatthew G. Knepley w[i*order+j] = 0.5 * wx[i] * wy[j]; 758494e7359SMatthew G. Knepley } 759494e7359SMatthew G. Knepley } 760494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 761494e7359SMatthew G. Knepley break; 762494e7359SMatthew G. Knepley case 3: 763dcca6d9dSJed Brown ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr); 764552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 765552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 766552aa4f7SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 767552aa4f7SMatthew G. Knepley for (i = 0; i < order; ++i) { 768552aa4f7SMatthew G. Knepley for (j = 0; j < order; ++j) { 769552aa4f7SMatthew G. Knepley for (k = 0; k < order; ++k) { 770552aa4f7SMatthew 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); 771552aa4f7SMatthew G. Knepley w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k]; 772494e7359SMatthew G. Knepley } 773494e7359SMatthew G. Knepley } 774494e7359SMatthew G. Knepley } 775494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 776494e7359SMatthew G. Knepley break; 777494e7359SMatthew G. Knepley default: 778494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 779494e7359SMatthew G. Knepley } 78021454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 781bcede257SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*q, order);CHKERRQ(ierr); 78221454ff5SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr); 783494e7359SMatthew G. Knepley PetscFunctionReturn(0); 784494e7359SMatthew G. Knepley } 785494e7359SMatthew G. Knepley 786494e7359SMatthew G. Knepley #undef __FUNCT__ 787b3c0f97bSTom Klotz #define __FUNCT__ "PetscDTTanhSinhTensorQuadrature" 788b3c0f97bSTom Klotz /*@C 789b3c0f97bSTom Klotz PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 790b3c0f97bSTom Klotz 791b3c0f97bSTom Klotz Not Collective 792b3c0f97bSTom Klotz 793b3c0f97bSTom Klotz Input Arguments: 794b3c0f97bSTom Klotz + dim - The cell dimension 795b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l 796b3c0f97bSTom Klotz . a - left end of interval (often-1) 797b3c0f97bSTom Klotz - b - right end of interval (often +1) 798b3c0f97bSTom Klotz 799b3c0f97bSTom Klotz Output Argument: 800b3c0f97bSTom Klotz . q - A PetscQuadrature object 801b3c0f97bSTom Klotz 802b3c0f97bSTom Klotz Level: intermediate 803b3c0f97bSTom Klotz 804b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature() 805b3c0f97bSTom Klotz @*/ 806b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 807b3c0f97bSTom Klotz { 808b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 809b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 810b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 811b3c0f97bSTom Klotz const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 812d84b4d08SMatthew G. Knepley PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 813b3c0f97bSTom Klotz PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 814b3c0f97bSTom Klotz PetscReal *x, *w; 815b3c0f97bSTom Klotz PetscInt K, k, npoints; 816b3c0f97bSTom Klotz PetscErrorCode ierr; 817b3c0f97bSTom Klotz 818b3c0f97bSTom Klotz PetscFunctionBegin; 819b3c0f97bSTom Klotz if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 820b3c0f97bSTom Klotz if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 821b3c0f97bSTom Klotz /* Find K such that the weights are < 32 digits of precision */ 822b3c0f97bSTom Klotz for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 8239add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 824b3c0f97bSTom Klotz } 825b3c0f97bSTom Klotz ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 826b3c0f97bSTom Klotz ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 827b3c0f97bSTom Klotz npoints = 2*K-1; 828b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 829b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 830b3c0f97bSTom Klotz /* Center term */ 831b3c0f97bSTom Klotz x[0] = beta; 832b3c0f97bSTom Klotz w[0] = 0.5*alpha*PETSC_PI; 833b3c0f97bSTom Klotz for (k = 1; k < K; ++k) { 8349add2064SThomas Klotz wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 8359add2064SThomas Klotz xk = tanh(0.5*PETSC_PI*PetscSinhReal(k*h)); 836b3c0f97bSTom Klotz x[2*k-1] = -alpha*xk+beta; 837b3c0f97bSTom Klotz w[2*k-1] = wk; 838b3c0f97bSTom Klotz x[2*k+0] = alpha*xk+beta; 839b3c0f97bSTom Klotz w[2*k+0] = wk; 840b3c0f97bSTom Klotz } 841b3c0f97bSTom Klotz ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr); 842b3c0f97bSTom Klotz PetscFunctionReturn(0); 843b3c0f97bSTom Klotz } 844b3c0f97bSTom Klotz 845b3c0f97bSTom Klotz #undef __FUNCT__ 846b3c0f97bSTom Klotz #define __FUNCT__ "PetscDTTanhSinhIntegrate" 847b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 848b3c0f97bSTom Klotz { 849b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 850b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 851b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 852b3c0f97bSTom Klotz PetscReal h = 1.0; /* Step size, length between x_k */ 853b3c0f97bSTom Klotz PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 854b3c0f97bSTom Klotz PetscReal osum = 0.0; /* Integral on last level */ 855b3c0f97bSTom Klotz PetscReal psum = 0.0; /* Integral on the level before the last level */ 856b3c0f97bSTom Klotz PetscReal sum; /* Integral on current level */ 857446c295cSMatthew G. Knepley PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 858b3c0f97bSTom Klotz PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 859b3c0f97bSTom Klotz PetscReal wk; /* Quadrature weight at x_k */ 860b3c0f97bSTom Klotz PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 861b3c0f97bSTom Klotz PetscInt d; /* Digits of precision in the integral */ 862b3c0f97bSTom Klotz 863b3c0f97bSTom Klotz PetscFunctionBegin; 864b3c0f97bSTom Klotz if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 865b3c0f97bSTom Klotz /* Center term */ 866b3c0f97bSTom Klotz func(beta, &lval); 867b3c0f97bSTom Klotz sum = 0.5*alpha*PETSC_PI*lval; 868b3c0f97bSTom Klotz /* */ 869b3c0f97bSTom Klotz do { 870b3c0f97bSTom Klotz PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 871b3c0f97bSTom Klotz PetscInt k = 1; 872b3c0f97bSTom Klotz 873b3c0f97bSTom Klotz ++l; 874b3c0f97bSTom Klotz /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 875b3c0f97bSTom Klotz /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 876b3c0f97bSTom Klotz psum = osum; 877b3c0f97bSTom Klotz osum = sum; 878b3c0f97bSTom Klotz h *= 0.5; 879b3c0f97bSTom Klotz sum *= 0.5; 880b3c0f97bSTom Klotz do { 8819add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 882446c295cSMatthew G. Knepley yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 883446c295cSMatthew G. Knepley lx = -alpha*(1.0 - yk)+beta; 884446c295cSMatthew G. Knepley rx = alpha*(1.0 - yk)+beta; 885b3c0f97bSTom Klotz func(lx, &lval); 886b3c0f97bSTom Klotz func(rx, &rval); 887b3c0f97bSTom Klotz lterm = alpha*wk*lval; 888b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 889b3c0f97bSTom Klotz sum += lterm; 890b3c0f97bSTom Klotz rterm = alpha*wk*rval; 891b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 892b3c0f97bSTom Klotz sum += rterm; 893b3c0f97bSTom Klotz ++k; 894b3c0f97bSTom Klotz /* Only need to evaluate every other point on refined levels */ 895b3c0f97bSTom Klotz if (l != 1) ++k; 8969add2064SThomas Klotz } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 897b3c0f97bSTom Klotz 898b3c0f97bSTom Klotz d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 899b3c0f97bSTom Klotz d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 900b3c0f97bSTom Klotz d3 = PetscLog10Real(maxTerm) - p; 90109d48545SBarry Smith if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 90209d48545SBarry Smith else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 903b3c0f97bSTom Klotz d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 9049add2064SThomas Klotz } while (d < digits && l < 12); 905b3c0f97bSTom Klotz *sol = sum; 906e510cb1fSThomas Klotz 907b3c0f97bSTom Klotz PetscFunctionReturn(0); 908b3c0f97bSTom Klotz } 909b3c0f97bSTom Klotz 91029f144ccSMatthew G. Knepley #ifdef PETSC_HAVE_MPFR 91129f144ccSMatthew G. Knepley #undef __FUNCT__ 91229f144ccSMatthew G. Knepley #define __FUNCT__ "PetscDTTanhSinhIntegrateMPFR" 91329f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 91429f144ccSMatthew G. Knepley { 915e510cb1fSThomas Klotz const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 91629f144ccSMatthew G. Knepley PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 91729f144ccSMatthew G. Knepley mpfr_t alpha; /* Half-width of the integration interval */ 91829f144ccSMatthew G. Knepley mpfr_t beta; /* Center of the integration interval */ 91929f144ccSMatthew G. Knepley mpfr_t h; /* Step size, length between x_k */ 92029f144ccSMatthew G. Knepley mpfr_t osum; /* Integral on last level */ 92129f144ccSMatthew G. Knepley mpfr_t psum; /* Integral on the level before the last level */ 92229f144ccSMatthew G. Knepley mpfr_t sum; /* Integral on current level */ 92329f144ccSMatthew G. Knepley mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 92429f144ccSMatthew G. Knepley mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 92529f144ccSMatthew G. Knepley mpfr_t wk; /* Quadrature weight at x_k */ 92629f144ccSMatthew G. Knepley PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 92729f144ccSMatthew G. Knepley PetscInt d; /* Digits of precision in the integral */ 92829f144ccSMatthew G. Knepley mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 92929f144ccSMatthew G. Knepley 93029f144ccSMatthew G. Knepley PetscFunctionBegin; 93129f144ccSMatthew G. Knepley if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 93229f144ccSMatthew G. Knepley /* Create high precision storage */ 933c9f744b5SMatthew G. Knepley mpfr_inits2(PetscCeilReal(safetyFactor*digits*PetscLogReal(10.)/PetscLogReal(2.)), alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 93429f144ccSMatthew G. Knepley /* Initialization */ 93529f144ccSMatthew G. Knepley mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 93629f144ccSMatthew G. Knepley mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 93729f144ccSMatthew G. Knepley mpfr_set_d(osum, 0.0, MPFR_RNDN); 93829f144ccSMatthew G. Knepley mpfr_set_d(psum, 0.0, MPFR_RNDN); 93929f144ccSMatthew G. Knepley mpfr_set_d(h, 1.0, MPFR_RNDN); 94029f144ccSMatthew G. Knepley mpfr_const_pi(pi2, MPFR_RNDN); 94129f144ccSMatthew G. Knepley mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 94229f144ccSMatthew G. Knepley /* Center term */ 94329f144ccSMatthew G. Knepley func(0.5*(b+a), &lval); 94429f144ccSMatthew G. Knepley mpfr_set(sum, pi2, MPFR_RNDN); 94529f144ccSMatthew G. Knepley mpfr_mul(sum, sum, alpha, MPFR_RNDN); 94629f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 94729f144ccSMatthew G. Knepley /* */ 94829f144ccSMatthew G. Knepley do { 94929f144ccSMatthew G. Knepley PetscReal d1, d2, d3, d4; 95029f144ccSMatthew G. Knepley PetscInt k = 1; 95129f144ccSMatthew G. Knepley 95229f144ccSMatthew G. Knepley ++l; 95329f144ccSMatthew G. Knepley mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 95429f144ccSMatthew G. Knepley /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 95529f144ccSMatthew G. Knepley /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 95629f144ccSMatthew G. Knepley mpfr_set(psum, osum, MPFR_RNDN); 95729f144ccSMatthew G. Knepley mpfr_set(osum, sum, MPFR_RNDN); 95829f144ccSMatthew G. Knepley mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 95929f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 96029f144ccSMatthew G. Knepley do { 96129f144ccSMatthew G. Knepley mpfr_set_si(kh, k, MPFR_RNDN); 96229f144ccSMatthew G. Knepley mpfr_mul(kh, kh, h, MPFR_RNDN); 96329f144ccSMatthew G. Knepley /* Weight */ 96429f144ccSMatthew G. Knepley mpfr_set(wk, h, MPFR_RNDN); 96529f144ccSMatthew G. Knepley mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 96629f144ccSMatthew G. Knepley mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 96729f144ccSMatthew G. Knepley mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 96829f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 96929f144ccSMatthew G. Knepley mpfr_sqr(tmp, tmp, MPFR_RNDN); 97029f144ccSMatthew G. Knepley mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 97129f144ccSMatthew G. Knepley mpfr_div(wk, wk, tmp, MPFR_RNDN); 97229f144ccSMatthew G. Knepley /* Abscissa */ 97329f144ccSMatthew G. Knepley mpfr_set_d(yk, 1.0, MPFR_RNDZ); 97429f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 97529f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 97629f144ccSMatthew G. Knepley mpfr_exp(tmp, msinh, MPFR_RNDN); 97729f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 97829f144ccSMatthew G. Knepley /* Quadrature points */ 97929f144ccSMatthew G. Knepley mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 98029f144ccSMatthew G. Knepley mpfr_mul(lx, lx, alpha, MPFR_RNDU); 98129f144ccSMatthew G. Knepley mpfr_add(lx, lx, beta, MPFR_RNDU); 98229f144ccSMatthew G. Knepley mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 98329f144ccSMatthew G. Knepley mpfr_mul(rx, rx, alpha, MPFR_RNDD); 98429f144ccSMatthew G. Knepley mpfr_add(rx, rx, beta, MPFR_RNDD); 98529f144ccSMatthew G. Knepley /* Evaluation */ 98629f144ccSMatthew G. Knepley func(mpfr_get_d(lx, MPFR_RNDU), &lval); 98729f144ccSMatthew G. Knepley func(mpfr_get_d(rx, MPFR_RNDD), &rval); 98829f144ccSMatthew G. Knepley /* Update */ 98929f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 99029f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 99129f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 99229f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 99329f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 99429f144ccSMatthew G. Knepley mpfr_set(curTerm, tmp, MPFR_RNDN); 99529f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 99629f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 99729f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 99829f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 99929f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 100029f144ccSMatthew G. Knepley mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 100129f144ccSMatthew G. Knepley ++k; 100229f144ccSMatthew G. Knepley /* Only need to evaluate every other point on refined levels */ 100329f144ccSMatthew G. Knepley if (l != 1) ++k; 100429f144ccSMatthew G. Knepley mpfr_log10(tmp, wk, MPFR_RNDN); 100529f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 1006c9f744b5SMatthew G. Knepley } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 100729f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, osum, MPFR_RNDN); 100829f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 100929f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 101029f144ccSMatthew G. Knepley d1 = mpfr_get_d(tmp, MPFR_RNDN); 101129f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, psum, MPFR_RNDN); 101229f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 101329f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 101429f144ccSMatthew G. Knepley d2 = mpfr_get_d(tmp, MPFR_RNDN); 101529f144ccSMatthew G. Knepley mpfr_log10(tmp, maxTerm, MPFR_RNDN); 1016c9f744b5SMatthew G. Knepley d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 101729f144ccSMatthew G. Knepley mpfr_log10(tmp, curTerm, MPFR_RNDN); 101829f144ccSMatthew G. Knepley d4 = mpfr_get_d(tmp, MPFR_RNDN); 101929f144ccSMatthew G. Knepley d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1020b0649871SThomas Klotz } while (d < digits && l < 8); 102129f144ccSMatthew G. Knepley *sol = mpfr_get_d(sum, MPFR_RNDN); 102229f144ccSMatthew G. Knepley /* Cleanup */ 102329f144ccSMatthew G. Knepley mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 102429f144ccSMatthew G. Knepley PetscFunctionReturn(0); 102529f144ccSMatthew G. Knepley } 1026d525116cSMatthew G. Knepley #else 1027d525116cSMatthew G. Knepley #undef __FUNCT__ 1028d525116cSMatthew G. Knepley #define __FUNCT__ "PetscDTTanhSinhIntegrateMPFR" 1029d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1030d525116cSMatthew G. Knepley { 1031d525116cSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 1032d525116cSMatthew G. Knepley } 103329f144ccSMatthew G. Knepley #endif 103429f144ccSMatthew G. Knepley 1035b3c0f97bSTom Klotz #undef __FUNCT__ 1036194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR" 1037194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 1038194825f6SJed Brown * A in column-major format 1039194825f6SJed Brown * Ainv in row-major format 1040194825f6SJed Brown * tau has length m 1041194825f6SJed Brown * worksize must be >= max(1,n) 1042194825f6SJed Brown */ 1043194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1044194825f6SJed Brown { 1045194825f6SJed Brown PetscErrorCode ierr; 1046194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1047194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 1048194825f6SJed Brown 1049194825f6SJed Brown PetscFunctionBegin; 1050194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1051194825f6SJed Brown { 1052194825f6SJed Brown PetscInt i,j; 1053dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 1054194825f6SJed Brown for (j=0; j<n; j++) { 1055194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 1056194825f6SJed Brown } 1057194825f6SJed Brown mstride = m; 1058194825f6SJed Brown } 1059194825f6SJed Brown #else 1060194825f6SJed Brown A = A_in; 1061194825f6SJed Brown Ainv = Ainv_out; 1062194825f6SJed Brown #endif 1063194825f6SJed Brown 1064194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 1065194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 1066194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 1067194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 1068194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1069001a771dSBarry Smith PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1070194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 1071194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1072194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1073194825f6SJed Brown 1074194825f6SJed Brown /* Extract an explicit representation of Q */ 1075194825f6SJed Brown Q = Ainv; 1076194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 1077194825f6SJed Brown K = N; /* full rank */ 1078001a771dSBarry Smith PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1079194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1080194825f6SJed Brown 1081194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1082194825f6SJed Brown Alpha = 1.0; 1083194825f6SJed Brown ldb = lda; 1084001a771dSBarry Smith PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 1085194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 1086194825f6SJed Brown 1087194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1088194825f6SJed Brown { 1089194825f6SJed Brown PetscInt i; 1090194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 1091194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 1092194825f6SJed Brown } 1093194825f6SJed Brown #endif 1094194825f6SJed Brown PetscFunctionReturn(0); 1095194825f6SJed Brown } 1096194825f6SJed Brown 1097194825f6SJed Brown #undef __FUNCT__ 1098194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate" 1099194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 1100194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 1101194825f6SJed Brown { 1102194825f6SJed Brown PetscErrorCode ierr; 1103194825f6SJed Brown PetscReal *Bv; 1104194825f6SJed Brown PetscInt i,j; 1105194825f6SJed Brown 1106194825f6SJed Brown PetscFunctionBegin; 1107785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 1108194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 1109194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 1110194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 1111194825f6SJed Brown for (i=0; i<ninterval; i++) { 1112194825f6SJed Brown for (j=0; j<ndegree; j++) { 1113194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1114194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1115194825f6SJed Brown } 1116194825f6SJed Brown } 1117194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 1118194825f6SJed Brown PetscFunctionReturn(0); 1119194825f6SJed Brown } 1120194825f6SJed Brown 1121194825f6SJed Brown #undef __FUNCT__ 1122194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly" 1123194825f6SJed Brown /*@ 1124194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 1125194825f6SJed Brown 1126194825f6SJed Brown Not Collective 1127194825f6SJed Brown 1128194825f6SJed Brown Input Arguments: 1129194825f6SJed Brown + degree - degree of reconstruction polynomial 1130194825f6SJed Brown . nsource - number of source intervals 1131194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 1132194825f6SJed Brown . ntarget - number of target intervals 1133194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 1134194825f6SJed Brown 1135194825f6SJed Brown Output Arguments: 1136194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 1137194825f6SJed Brown 1138194825f6SJed Brown Level: advanced 1139194825f6SJed Brown 1140194825f6SJed Brown .seealso: PetscDTLegendreEval() 1141194825f6SJed Brown @*/ 1142194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 1143194825f6SJed Brown { 1144194825f6SJed Brown PetscErrorCode ierr; 1145194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 1146194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 1147194825f6SJed Brown PetscScalar *tau,*work; 1148194825f6SJed Brown 1149194825f6SJed Brown PetscFunctionBegin; 1150194825f6SJed Brown PetscValidRealPointer(sourcex,3); 1151194825f6SJed Brown PetscValidRealPointer(targetx,5); 1152194825f6SJed Brown PetscValidRealPointer(R,6); 1153194825f6SJed 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); 1154194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 1155194825f6SJed Brown for (i=0; i<nsource; i++) { 115657622a8eSBarry 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]); 1157194825f6SJed Brown } 1158194825f6SJed Brown for (i=0; i<ntarget; i++) { 115957622a8eSBarry 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]); 1160194825f6SJed Brown } 1161194825f6SJed Brown #endif 1162194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 1163194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 1164194825f6SJed Brown center = (xmin + xmax)/2; 1165194825f6SJed Brown hscale = (xmax - xmin)/2; 1166194825f6SJed Brown worksize = nsource; 1167dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 1168dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 1169194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 1170194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 1171194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 1172194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 1173194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 1174194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 1175194825f6SJed Brown for (i=0; i<ntarget; i++) { 1176194825f6SJed Brown PetscReal rowsum = 0; 1177194825f6SJed Brown for (j=0; j<nsource; j++) { 1178194825f6SJed Brown PetscReal sum = 0; 1179194825f6SJed Brown for (k=0; k<degree+1; k++) { 1180194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 1181194825f6SJed Brown } 1182194825f6SJed Brown R[i*nsource+j] = sum; 1183194825f6SJed Brown rowsum += sum; 1184194825f6SJed Brown } 1185194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 1186194825f6SJed Brown } 1187194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 1188194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 1189194825f6SJed Brown PetscFunctionReturn(0); 1190194825f6SJed Brown } 1191