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 7*497880caSRichard Tran Mills #if defined(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 2940d8ff71SMatthew G. Knepley /*@ 3040d8ff71SMatthew G. Knepley PetscQuadratureCreate - Create a PetscQuadrature object 3140d8ff71SMatthew G. Knepley 3240d8ff71SMatthew G. Knepley Collective on MPI_Comm 3340d8ff71SMatthew G. Knepley 3440d8ff71SMatthew G. Knepley Input Parameter: 3540d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object 3640d8ff71SMatthew G. Knepley 3740d8ff71SMatthew G. Knepley Output Parameter: 3840d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 3940d8ff71SMatthew G. Knepley 4040d8ff71SMatthew G. Knepley Level: beginner 4140d8ff71SMatthew G. Knepley 4240d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, create 4340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 4440d8ff71SMatthew G. Knepley @*/ 4521454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 4621454ff5SMatthew G. Knepley { 4721454ff5SMatthew G. Knepley PetscErrorCode ierr; 4821454ff5SMatthew G. Knepley 4921454ff5SMatthew G. Knepley PetscFunctionBegin; 5021454ff5SMatthew G. Knepley PetscValidPointer(q, 2); 51623436dcSMatthew G. Knepley ierr = PetscSysInitializePackage();CHKERRQ(ierr); 5273107ff1SLisandro Dalcin ierr = PetscHeaderCreate(*q,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 5321454ff5SMatthew G. Knepley (*q)->dim = -1; 54a6b92713SMatthew G. Knepley (*q)->Nc = 1; 55bcede257SMatthew G. Knepley (*q)->order = -1; 5621454ff5SMatthew G. Knepley (*q)->numPoints = 0; 5721454ff5SMatthew G. Knepley (*q)->points = NULL; 5821454ff5SMatthew G. Knepley (*q)->weights = NULL; 5921454ff5SMatthew G. Knepley PetscFunctionReturn(0); 6021454ff5SMatthew G. Knepley } 6121454ff5SMatthew G. Knepley 62c9638911SMatthew G. Knepley /*@ 63c9638911SMatthew G. Knepley PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 64c9638911SMatthew G. Knepley 65c9638911SMatthew G. Knepley Collective on PetscQuadrature 66c9638911SMatthew G. Knepley 67c9638911SMatthew G. Knepley Input Parameter: 68c9638911SMatthew G. Knepley . q - The PetscQuadrature object 69c9638911SMatthew G. Knepley 70c9638911SMatthew G. Knepley Output Parameter: 71c9638911SMatthew G. Knepley . r - The new PetscQuadrature object 72c9638911SMatthew G. Knepley 73c9638911SMatthew G. Knepley Level: beginner 74c9638911SMatthew G. Knepley 75c9638911SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, clone 76c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 77c9638911SMatthew G. Knepley @*/ 78c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 79c9638911SMatthew G. Knepley { 80a6b92713SMatthew G. Knepley PetscInt order, dim, Nc, Nq; 81c9638911SMatthew G. Knepley const PetscReal *points, *weights; 82c9638911SMatthew G. Knepley PetscReal *p, *w; 83c9638911SMatthew G. Knepley PetscErrorCode ierr; 84c9638911SMatthew G. Knepley 85c9638911SMatthew G. Knepley PetscFunctionBegin; 86c9638911SMatthew G. Knepley PetscValidPointer(q, 2); 87c9638911SMatthew G. Knepley ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 88c9638911SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 89c9638911SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 90a6b92713SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 91c9638911SMatthew G. Knepley ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 92f0a0bfafSMatthew G. Knepley ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr); 93c9638911SMatthew G. Knepley ierr = PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));CHKERRQ(ierr); 94a6b92713SMatthew G. Knepley ierr = PetscMemcpy(w, weights, Nc * Nq * sizeof(PetscReal));CHKERRQ(ierr); 95a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr); 96c9638911SMatthew G. Knepley PetscFunctionReturn(0); 97c9638911SMatthew G. Knepley } 98c9638911SMatthew G. Knepley 9940d8ff71SMatthew G. Knepley /*@ 10040d8ff71SMatthew G. Knepley PetscQuadratureDestroy - Destroys a PetscQuadrature object 10140d8ff71SMatthew G. Knepley 10240d8ff71SMatthew G. Knepley Collective on PetscQuadrature 10340d8ff71SMatthew G. Knepley 10440d8ff71SMatthew G. Knepley Input Parameter: 10540d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 10640d8ff71SMatthew G. Knepley 10740d8ff71SMatthew G. Knepley Level: beginner 10840d8ff71SMatthew G. Knepley 10940d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, destroy 11040d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 11140d8ff71SMatthew G. Knepley @*/ 112bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 113bfa639d9SMatthew G. Knepley { 114bfa639d9SMatthew G. Knepley PetscErrorCode ierr; 115bfa639d9SMatthew G. Knepley 116bfa639d9SMatthew G. Knepley PetscFunctionBegin; 11721454ff5SMatthew G. Knepley if (!*q) PetscFunctionReturn(0); 11821454ff5SMatthew G. Knepley PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 11921454ff5SMatthew G. Knepley if (--((PetscObject)(*q))->refct > 0) { 12021454ff5SMatthew G. Knepley *q = NULL; 12121454ff5SMatthew G. Knepley PetscFunctionReturn(0); 12221454ff5SMatthew G. Knepley } 12321454ff5SMatthew G. Knepley ierr = PetscFree((*q)->points);CHKERRQ(ierr); 12421454ff5SMatthew G. Knepley ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 12521454ff5SMatthew G. Knepley ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 12621454ff5SMatthew G. Knepley PetscFunctionReturn(0); 12721454ff5SMatthew G. Knepley } 12821454ff5SMatthew G. Knepley 129bcede257SMatthew G. Knepley /*@ 130a6b92713SMatthew G. Knepley PetscQuadratureGetOrder - Return the order of the method 131bcede257SMatthew G. Knepley 132bcede257SMatthew G. Knepley Not collective 133bcede257SMatthew G. Knepley 134bcede257SMatthew G. Knepley Input Parameter: 135bcede257SMatthew G. Knepley . q - The PetscQuadrature object 136bcede257SMatthew G. Knepley 137bcede257SMatthew G. Knepley Output Parameter: 138bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 139bcede257SMatthew G. Knepley 140bcede257SMatthew G. Knepley Level: intermediate 141bcede257SMatthew G. Knepley 142bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 143bcede257SMatthew G. Knepley @*/ 144bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 145bcede257SMatthew G. Knepley { 146bcede257SMatthew G. Knepley PetscFunctionBegin; 147bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 148bcede257SMatthew G. Knepley PetscValidPointer(order, 2); 149bcede257SMatthew G. Knepley *order = q->order; 150bcede257SMatthew G. Knepley PetscFunctionReturn(0); 151bcede257SMatthew G. Knepley } 152bcede257SMatthew G. Knepley 153bcede257SMatthew G. Knepley /*@ 154a6b92713SMatthew G. Knepley PetscQuadratureSetOrder - Return the order of the method 155bcede257SMatthew G. Knepley 156bcede257SMatthew G. Knepley Not collective 157bcede257SMatthew G. Knepley 158bcede257SMatthew G. Knepley Input Parameters: 159bcede257SMatthew G. Knepley + q - The PetscQuadrature object 160bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 161bcede257SMatthew G. Knepley 162bcede257SMatthew G. Knepley Level: intermediate 163bcede257SMatthew G. Knepley 164bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 165bcede257SMatthew G. Knepley @*/ 166bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 167bcede257SMatthew G. Knepley { 168bcede257SMatthew G. Knepley PetscFunctionBegin; 169bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 170bcede257SMatthew G. Knepley q->order = order; 171bcede257SMatthew G. Knepley PetscFunctionReturn(0); 172bcede257SMatthew G. Knepley } 173bcede257SMatthew G. Knepley 174a6b92713SMatthew G. Knepley /*@ 175a6b92713SMatthew G. Knepley PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated 176a6b92713SMatthew G. Knepley 177a6b92713SMatthew G. Knepley Not collective 178a6b92713SMatthew G. Knepley 179a6b92713SMatthew G. Knepley Input Parameter: 180a6b92713SMatthew G. Knepley . q - The PetscQuadrature object 181a6b92713SMatthew G. Knepley 182a6b92713SMatthew G. Knepley Output Parameter: 183a6b92713SMatthew G. Knepley . Nc - The number of components 184a6b92713SMatthew G. Knepley 185a6b92713SMatthew G. Knepley Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 186a6b92713SMatthew G. Knepley 187a6b92713SMatthew G. Knepley Level: intermediate 188a6b92713SMatthew G. Knepley 189a6b92713SMatthew G. Knepley .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 190a6b92713SMatthew G. Knepley @*/ 191a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc) 192a6b92713SMatthew G. Knepley { 193a6b92713SMatthew G. Knepley PetscFunctionBegin; 194a6b92713SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 195a6b92713SMatthew G. Knepley PetscValidPointer(Nc, 2); 196a6b92713SMatthew G. Knepley *Nc = q->Nc; 197a6b92713SMatthew G. Knepley PetscFunctionReturn(0); 198a6b92713SMatthew G. Knepley } 199a6b92713SMatthew G. Knepley 200a6b92713SMatthew G. Knepley /*@ 201a6b92713SMatthew G. Knepley PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated 202a6b92713SMatthew G. Knepley 203a6b92713SMatthew G. Knepley Not collective 204a6b92713SMatthew G. Knepley 205a6b92713SMatthew G. Knepley Input Parameters: 206a6b92713SMatthew G. Knepley + q - The PetscQuadrature object 207a6b92713SMatthew G. Knepley - Nc - The number of components 208a6b92713SMatthew G. Knepley 209a6b92713SMatthew G. Knepley Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 210a6b92713SMatthew G. Knepley 211a6b92713SMatthew G. Knepley Level: intermediate 212a6b92713SMatthew G. Knepley 213a6b92713SMatthew G. Knepley .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 214a6b92713SMatthew G. Knepley @*/ 215a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc) 216a6b92713SMatthew G. Knepley { 217a6b92713SMatthew G. Knepley PetscFunctionBegin; 218a6b92713SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 219a6b92713SMatthew G. Knepley q->Nc = Nc; 220a6b92713SMatthew G. Knepley PetscFunctionReturn(0); 221a6b92713SMatthew G. Knepley } 222a6b92713SMatthew G. Knepley 22340d8ff71SMatthew G. Knepley /*@C 22440d8ff71SMatthew G. Knepley PetscQuadratureGetData - Returns the data defining the quadrature 22540d8ff71SMatthew G. Knepley 22640d8ff71SMatthew G. Knepley Not collective 22740d8ff71SMatthew G. Knepley 22840d8ff71SMatthew G. Knepley Input Parameter: 22940d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 23040d8ff71SMatthew G. Knepley 23140d8ff71SMatthew G. Knepley Output Parameters: 23240d8ff71SMatthew G. Knepley + dim - The spatial dimension 233805e7170SToby Isaac . Nc - The number of components 23440d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 23540d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 23640d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 23740d8ff71SMatthew G. Knepley 23840d8ff71SMatthew G. Knepley Level: intermediate 23940d8ff71SMatthew G. Knepley 24095452b02SPatrick Sanan Fortran Notes: 24195452b02SPatrick Sanan From Fortran you must call PetscQuadratureRestoreData() when you are done with the data 2421fd49c25SBarry Smith 24340d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature 24440d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 24540d8ff71SMatthew G. Knepley @*/ 246a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 24721454ff5SMatthew G. Knepley { 24821454ff5SMatthew G. Knepley PetscFunctionBegin; 24921454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 25021454ff5SMatthew G. Knepley if (dim) { 25121454ff5SMatthew G. Knepley PetscValidPointer(dim, 2); 25221454ff5SMatthew G. Knepley *dim = q->dim; 25321454ff5SMatthew G. Knepley } 254a6b92713SMatthew G. Knepley if (Nc) { 255a6b92713SMatthew G. Knepley PetscValidPointer(Nc, 3); 256a6b92713SMatthew G. Knepley *Nc = q->Nc; 257a6b92713SMatthew G. Knepley } 25821454ff5SMatthew G. Knepley if (npoints) { 259a6b92713SMatthew G. Knepley PetscValidPointer(npoints, 4); 26021454ff5SMatthew G. Knepley *npoints = q->numPoints; 26121454ff5SMatthew G. Knepley } 26221454ff5SMatthew G. Knepley if (points) { 263a6b92713SMatthew G. Knepley PetscValidPointer(points, 5); 26421454ff5SMatthew G. Knepley *points = q->points; 26521454ff5SMatthew G. Knepley } 26621454ff5SMatthew G. Knepley if (weights) { 267a6b92713SMatthew G. Knepley PetscValidPointer(weights, 6); 26821454ff5SMatthew G. Knepley *weights = q->weights; 26921454ff5SMatthew G. Knepley } 27021454ff5SMatthew G. Knepley PetscFunctionReturn(0); 27121454ff5SMatthew G. Knepley } 27221454ff5SMatthew G. Knepley 27340d8ff71SMatthew G. Knepley /*@C 27440d8ff71SMatthew G. Knepley PetscQuadratureSetData - Sets the data defining the quadrature 27540d8ff71SMatthew G. Knepley 27640d8ff71SMatthew G. Knepley Not collective 27740d8ff71SMatthew G. Knepley 27840d8ff71SMatthew G. Knepley Input Parameters: 27940d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 28040d8ff71SMatthew G. Knepley . dim - The spatial dimension 281a6b92713SMatthew G. Knepley , Nc - The number of components 28240d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 28340d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 28440d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 28540d8ff71SMatthew G. Knepley 286f2fd9e53SMatthew G. Knepley Note: This routine owns the references to points and weights, so they msut be allocated using PetscMalloc() and the user should not free them. 287f2fd9e53SMatthew G. Knepley 28840d8ff71SMatthew G. Knepley Level: intermediate 28940d8ff71SMatthew G. Knepley 29040d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature 29140d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 29240d8ff71SMatthew G. Knepley @*/ 293a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 29421454ff5SMatthew G. Knepley { 29521454ff5SMatthew G. Knepley PetscFunctionBegin; 29621454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 29721454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 298a6b92713SMatthew G. Knepley if (Nc >= 0) q->Nc = Nc; 29921454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 30021454ff5SMatthew G. Knepley if (points) { 30121454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 30221454ff5SMatthew G. Knepley q->points = points; 30321454ff5SMatthew G. Knepley } 30421454ff5SMatthew G. Knepley if (weights) { 30521454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 30621454ff5SMatthew G. Knepley q->weights = weights; 30721454ff5SMatthew G. Knepley } 308f9fd7fdbSMatthew G. Knepley PetscFunctionReturn(0); 309f9fd7fdbSMatthew G. Knepley } 310f9fd7fdbSMatthew G. Knepley 31140d8ff71SMatthew G. Knepley /*@C 31240d8ff71SMatthew G. Knepley PetscQuadratureView - Views a PetscQuadrature object 31340d8ff71SMatthew G. Knepley 31440d8ff71SMatthew G. Knepley Collective on PetscQuadrature 31540d8ff71SMatthew G. Knepley 31640d8ff71SMatthew G. Knepley Input Parameters: 31740d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 31840d8ff71SMatthew G. Knepley - viewer - The PetscViewer object 31940d8ff71SMatthew G. Knepley 32040d8ff71SMatthew G. Knepley Level: beginner 32140d8ff71SMatthew G. Knepley 32240d8ff71SMatthew G. Knepley .keywords: PetscQuadrature, quadrature, view 32340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 32440d8ff71SMatthew G. Knepley @*/ 325f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 326f9fd7fdbSMatthew G. Knepley { 327a6b92713SMatthew G. Knepley PetscInt q, d, c; 328f9fd7fdbSMatthew G. Knepley PetscErrorCode ierr; 329f9fd7fdbSMatthew G. Knepley 330f9fd7fdbSMatthew G. Knepley PetscFunctionBegin; 33198c3331eSBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr); 332a6b92713SMatthew G. Knepley if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %D points with %D components\n (", quad->numPoints, quad->Nc);CHKERRQ(ierr);} 333a6b92713SMatthew G. Knepley else {ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %D points\n (", quad->numPoints);CHKERRQ(ierr);} 33421454ff5SMatthew G. Knepley for (q = 0; q < quad->numPoints; ++q) { 33521454ff5SMatthew G. Knepley for (d = 0; d < quad->dim; ++d) { 336f9fd7fdbSMatthew G. Knepley if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 337ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 338f9fd7fdbSMatthew G. Knepley } 339a6b92713SMatthew G. Knepley if (quad->Nc > 1) { 340a6b92713SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, ") (");CHKERRQ(ierr); 341a6b92713SMatthew G. Knepley for (c = 0; c < quad->Nc; ++c) { 342a6b92713SMatthew G. Knepley if (c) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 343a6b92713SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr); 344a6b92713SMatthew G. Knepley } 345a6b92713SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, ")\n");CHKERRQ(ierr); 346a6b92713SMatthew G. Knepley } else { 347ab15ae43SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr); 348f9fd7fdbSMatthew G. Knepley } 349a6b92713SMatthew G. Knepley } 350bfa639d9SMatthew G. Knepley PetscFunctionReturn(0); 351bfa639d9SMatthew G. Knepley } 352bfa639d9SMatthew G. Knepley 35389710940SMatthew G. Knepley /*@C 35489710940SMatthew G. Knepley PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 35589710940SMatthew G. Knepley 35689710940SMatthew G. Knepley Not collective 35789710940SMatthew G. Knepley 35889710940SMatthew G. Knepley Input Parameter: 35989710940SMatthew G. Knepley + q - The original PetscQuadrature 36089710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into 36189710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement 36289710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement 36389710940SMatthew G. Knepley 36489710940SMatthew G. Knepley Output Parameters: 36589710940SMatthew G. Knepley . dim - The dimension 36689710940SMatthew G. Knepley 36789710940SMatthew G. Knepley Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 36889710940SMatthew G. Knepley 369f5f57ec0SBarry Smith Not available from Fortran 370f5f57ec0SBarry Smith 37189710940SMatthew G. Knepley Level: intermediate 37289710940SMatthew G. Knepley 37389710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 37489710940SMatthew G. Knepley @*/ 37589710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 37689710940SMatthew G. Knepley { 37789710940SMatthew G. Knepley const PetscReal *points, *weights; 37889710940SMatthew G. Knepley PetscReal *pointsRef, *weightsRef; 379a6b92713SMatthew G. Knepley PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e; 38089710940SMatthew G. Knepley PetscErrorCode ierr; 38189710940SMatthew G. Knepley 38289710940SMatthew G. Knepley PetscFunctionBegin; 38389710940SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 38489710940SMatthew G. Knepley PetscValidPointer(v0, 3); 38589710940SMatthew G. Knepley PetscValidPointer(jac, 4); 38689710940SMatthew G. Knepley PetscValidPointer(qref, 5); 38789710940SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 38889710940SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 389a6b92713SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 39089710940SMatthew G. Knepley npointsRef = npoints*numSubelements; 39189710940SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 392a6b92713SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr); 39389710940SMatthew G. Knepley for (c = 0; c < numSubelements; ++c) { 39489710940SMatthew G. Knepley for (p = 0; p < npoints; ++p) { 39589710940SMatthew G. Knepley for (d = 0; d < dim; ++d) { 39689710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 39789710940SMatthew G. Knepley for (e = 0; e < dim; ++e) { 39889710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 39989710940SMatthew G. Knepley } 40089710940SMatthew G. Knepley } 40189710940SMatthew G. Knepley /* Could also use detJ here */ 402a6b92713SMatthew G. Knepley for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements; 40389710940SMatthew G. Knepley } 40489710940SMatthew G. Knepley } 40589710940SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 406a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 40789710940SMatthew G. Knepley PetscFunctionReturn(0); 40889710940SMatthew G. Knepley } 40989710940SMatthew G. Knepley 41037045ce4SJed Brown /*@ 41137045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 41237045ce4SJed Brown 41337045ce4SJed Brown Not Collective 41437045ce4SJed Brown 41537045ce4SJed Brown Input Arguments: 41637045ce4SJed Brown + npoints - number of spatial points to evaluate at 41737045ce4SJed Brown . points - array of locations to evaluate at 41837045ce4SJed Brown . ndegree - number of basis degrees to evaluate 41937045ce4SJed Brown - degrees - sorted array of degrees to evaluate 42037045ce4SJed Brown 42137045ce4SJed Brown Output Arguments: 4220298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 4230298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 4240298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 42537045ce4SJed Brown 42637045ce4SJed Brown Level: intermediate 42737045ce4SJed Brown 42837045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 42937045ce4SJed Brown @*/ 43037045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 43137045ce4SJed Brown { 43237045ce4SJed Brown PetscInt i,maxdegree; 43337045ce4SJed Brown 43437045ce4SJed Brown PetscFunctionBegin; 43537045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 43637045ce4SJed Brown maxdegree = degrees[ndegree-1]; 43737045ce4SJed Brown for (i=0; i<npoints; i++) { 43837045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 43937045ce4SJed Brown PetscInt j,k; 44037045ce4SJed Brown x = points[i]; 44137045ce4SJed Brown pm2 = 0; 44237045ce4SJed Brown pm1 = 1; 44337045ce4SJed Brown pd2 = 0; 44437045ce4SJed Brown pd1 = 0; 44537045ce4SJed Brown pdd2 = 0; 44637045ce4SJed Brown pdd1 = 0; 44737045ce4SJed Brown k = 0; 44837045ce4SJed Brown if (degrees[k] == 0) { 44937045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 45037045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 45137045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 45237045ce4SJed Brown k++; 45337045ce4SJed Brown } 45437045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 45537045ce4SJed Brown PetscReal p,d,dd; 45637045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 45737045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 45837045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 45937045ce4SJed Brown pm2 = pm1; 46037045ce4SJed Brown pm1 = p; 46137045ce4SJed Brown pd2 = pd1; 46237045ce4SJed Brown pd1 = d; 46337045ce4SJed Brown pdd2 = pdd1; 46437045ce4SJed Brown pdd1 = dd; 46537045ce4SJed Brown if (degrees[k] == j) { 46637045ce4SJed Brown if (B) B[i*ndegree+k] = p; 46737045ce4SJed Brown if (D) D[i*ndegree+k] = d; 46837045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 46937045ce4SJed Brown } 47037045ce4SJed Brown } 47137045ce4SJed Brown } 47237045ce4SJed Brown PetscFunctionReturn(0); 47337045ce4SJed Brown } 47437045ce4SJed Brown 47537045ce4SJed Brown /*@ 47637045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 47737045ce4SJed Brown 47837045ce4SJed Brown Not Collective 47937045ce4SJed Brown 48037045ce4SJed Brown Input Arguments: 48137045ce4SJed Brown + npoints - number of points 48237045ce4SJed Brown . a - left end of interval (often-1) 48337045ce4SJed Brown - b - right end of interval (often +1) 48437045ce4SJed Brown 48537045ce4SJed Brown Output Arguments: 48637045ce4SJed Brown + x - quadrature points 48737045ce4SJed Brown - w - quadrature weights 48837045ce4SJed Brown 48937045ce4SJed Brown Level: intermediate 49037045ce4SJed Brown 49137045ce4SJed Brown References: 49296a0c994SBarry Smith . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 49337045ce4SJed Brown 49437045ce4SJed Brown .seealso: PetscDTLegendreEval() 49537045ce4SJed Brown @*/ 49637045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 49737045ce4SJed Brown { 49837045ce4SJed Brown PetscErrorCode ierr; 49937045ce4SJed Brown PetscInt i; 50037045ce4SJed Brown PetscReal *work; 50137045ce4SJed Brown PetscScalar *Z; 50237045ce4SJed Brown PetscBLASInt N,LDZ,info; 50337045ce4SJed Brown 50437045ce4SJed Brown PetscFunctionBegin; 5050bfcf5a5SMatthew G. Knepley ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 50637045ce4SJed Brown /* Set up the Golub-Welsch system */ 50737045ce4SJed Brown for (i=0; i<npoints; i++) { 50837045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 50937045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 51037045ce4SJed Brown } 511dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 512c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 51337045ce4SJed Brown LDZ = N; 51437045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 5158b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 51637045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 5171c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 51837045ce4SJed Brown 51937045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 52037045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 52137045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 52219a57d60SBarry Smith if (x[i] == -0.0) x[i] = 0.0; 52337045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 5240d644c17SKarl Rupp 52588393a60SJed Brown w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 52637045ce4SJed Brown } 52737045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 52837045ce4SJed Brown PetscFunctionReturn(0); 52937045ce4SJed Brown } 530194825f6SJed Brown 531744bafbcSMatthew G. Knepley /*@ 532744bafbcSMatthew G. Knepley PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 533744bafbcSMatthew G. Knepley 534744bafbcSMatthew G. Knepley Not Collective 535744bafbcSMatthew G. Knepley 536744bafbcSMatthew G. Knepley Input Arguments: 537744bafbcSMatthew G. Knepley + dim - The spatial dimension 538a6b92713SMatthew G. Knepley . Nc - The number of components 539744bafbcSMatthew G. Knepley . npoints - number of points in one dimension 540744bafbcSMatthew G. Knepley . a - left end of interval (often-1) 541744bafbcSMatthew G. Knepley - b - right end of interval (often +1) 542744bafbcSMatthew G. Knepley 543744bafbcSMatthew G. Knepley Output Argument: 544744bafbcSMatthew G. Knepley . q - A PetscQuadrature object 545744bafbcSMatthew G. Knepley 546744bafbcSMatthew G. Knepley Level: intermediate 547744bafbcSMatthew G. Knepley 548744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 549744bafbcSMatthew G. Knepley @*/ 550a6b92713SMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 551744bafbcSMatthew G. Knepley { 552a6b92713SMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c; 553744bafbcSMatthew G. Knepley PetscReal *x, *w, *xw, *ww; 554744bafbcSMatthew G. Knepley PetscErrorCode ierr; 555744bafbcSMatthew G. Knepley 556744bafbcSMatthew G. Knepley PetscFunctionBegin; 557744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 558a6b92713SMatthew G. Knepley ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr); 559744bafbcSMatthew G. Knepley /* Set up the Golub-Welsch system */ 560744bafbcSMatthew G. Knepley switch (dim) { 561744bafbcSMatthew G. Knepley case 0: 562744bafbcSMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 563744bafbcSMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 564744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 565a6b92713SMatthew G. Knepley ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 566744bafbcSMatthew G. Knepley x[0] = 0.0; 567a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[c] = 1.0; 568744bafbcSMatthew G. Knepley break; 569744bafbcSMatthew G. Knepley case 1: 570a6b92713SMatthew G. Knepley ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr); 571a6b92713SMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr); 572a6b92713SMatthew G. Knepley for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i]; 573a6b92713SMatthew G. Knepley ierr = PetscFree(ww);CHKERRQ(ierr); 574744bafbcSMatthew G. Knepley break; 575744bafbcSMatthew G. Knepley case 2: 576744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 577744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 578744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 579744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 580744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+0] = xw[i]; 581744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+1] = xw[j]; 582a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j]; 583744bafbcSMatthew G. Knepley } 584744bafbcSMatthew G. Knepley } 585744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 586744bafbcSMatthew G. Knepley break; 587744bafbcSMatthew G. Knepley case 3: 588744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 589744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 590744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 591744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 592744bafbcSMatthew G. Knepley for (k = 0; k < npoints; ++k) { 593744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 594744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 595744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 596a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k]; 597744bafbcSMatthew G. Knepley } 598744bafbcSMatthew G. Knepley } 599744bafbcSMatthew G. Knepley } 600744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 601744bafbcSMatthew G. Knepley break; 602744bafbcSMatthew G. Knepley default: 603744bafbcSMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 604744bafbcSMatthew G. Knepley } 605744bafbcSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 6062f5fb066SToby Isaac ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 607a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 608744bafbcSMatthew G. Knepley PetscFunctionReturn(0); 609744bafbcSMatthew G. Knepley } 610744bafbcSMatthew G. Knepley 611494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 612494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 613494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 614494e7359SMatthew G. Knepley { 615494e7359SMatthew G. Knepley PetscReal f = 1.0; 616494e7359SMatthew G. Knepley PetscInt i; 617494e7359SMatthew G. Knepley 618494e7359SMatthew G. Knepley PetscFunctionBegin; 619494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 620494e7359SMatthew G. Knepley *factorial = f; 621494e7359SMatthew G. Knepley PetscFunctionReturn(0); 622494e7359SMatthew G. Knepley } 623494e7359SMatthew G. Knepley 624494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 625494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 626494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 627494e7359SMatthew G. Knepley { 628494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 629494e7359SMatthew G. Knepley PetscInt k; 630494e7359SMatthew G. Knepley 631494e7359SMatthew G. Knepley PetscFunctionBegin; 632494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 633494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 634494e7359SMatthew G. Knepley apb = a + b; 635494e7359SMatthew G. Knepley pn2 = 1.0; 636494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 637494e7359SMatthew G. Knepley *P = 0.0; 638494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 639494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 640494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 641494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 642494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 643494e7359SMatthew G. Knepley 644494e7359SMatthew G. Knepley a2 = a2 / a1; 645494e7359SMatthew G. Knepley a3 = a3 / a1; 646494e7359SMatthew G. Knepley a4 = a4 / a1; 647494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 648494e7359SMatthew G. Knepley pn2 = pn1; 649494e7359SMatthew G. Knepley pn1 = *P; 650494e7359SMatthew G. Knepley } 651494e7359SMatthew G. Knepley PetscFunctionReturn(0); 652494e7359SMatthew G. Knepley } 653494e7359SMatthew G. Knepley 654494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 655494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 656494e7359SMatthew G. Knepley { 657494e7359SMatthew G. Knepley PetscReal nP; 658494e7359SMatthew G. Knepley PetscErrorCode ierr; 659494e7359SMatthew G. Knepley 660494e7359SMatthew G. Knepley PetscFunctionBegin; 661494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 662494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 663494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 664494e7359SMatthew G. Knepley PetscFunctionReturn(0); 665494e7359SMatthew G. Knepley } 666494e7359SMatthew G. Knepley 667494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 668494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 669494e7359SMatthew G. Knepley { 670494e7359SMatthew G. Knepley PetscFunctionBegin; 671494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 672494e7359SMatthew G. Knepley *eta = y; 673494e7359SMatthew G. Knepley PetscFunctionReturn(0); 674494e7359SMatthew G. Knepley } 675494e7359SMatthew G. Knepley 676494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 677494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 678494e7359SMatthew G. Knepley { 679494e7359SMatthew G. Knepley PetscFunctionBegin; 680494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 681494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 682494e7359SMatthew G. Knepley *zeta = z; 683494e7359SMatthew G. Knepley PetscFunctionReturn(0); 684494e7359SMatthew G. Knepley } 685494e7359SMatthew G. Knepley 686494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 687494e7359SMatthew G. Knepley { 688494e7359SMatthew G. Knepley PetscInt maxIter = 100; 689494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 690a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 691494e7359SMatthew G. Knepley PetscInt k; 692494e7359SMatthew G. Knepley PetscErrorCode ierr; 693494e7359SMatthew G. Knepley 694494e7359SMatthew G. Knepley PetscFunctionBegin; 695a8291ba1SSatish Balay 6968b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 697a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 6980646a658SBarry Smith a2 = PetscTGamma(a + npoints + 1); 6990646a658SBarry Smith a3 = PetscTGamma(b + npoints + 1); 7000646a658SBarry Smith a4 = PetscTGamma(a + b + npoints + 1); 701a8291ba1SSatish Balay #else 70229bcbfd0SToby Isaac { 703d24bbb91SToby Isaac PetscInt ia, ib; 70429bcbfd0SToby Isaac 705d24bbb91SToby Isaac ia = (PetscInt) a; 706d24bbb91SToby Isaac ib = (PetscInt) b; 707d24bbb91SToby 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 */ 708d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr); 709d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr); 710d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr); 71129bcbfd0SToby Isaac } else { 712a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 71329bcbfd0SToby Isaac } 71429bcbfd0SToby Isaac } 715a8291ba1SSatish Balay #endif 716a8291ba1SSatish Balay 717494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 718494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 719494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 720494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 721494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 7228b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 723494e7359SMatthew G. Knepley PetscInt j; 724494e7359SMatthew G. Knepley 725494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 726494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 727494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 728494e7359SMatthew G. Knepley PetscInt i; 729494e7359SMatthew G. Knepley 730494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 731494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 732494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 733494e7359SMatthew G. Knepley delta = f / (fp - f * s); 734494e7359SMatthew G. Knepley r = r - delta; 73577b4d14cSPeter Brune if (PetscAbsReal(delta) < eps) break; 736494e7359SMatthew G. Knepley } 737494e7359SMatthew G. Knepley x[k] = r; 738494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 739494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 740494e7359SMatthew G. Knepley } 741494e7359SMatthew G. Knepley PetscFunctionReturn(0); 742494e7359SMatthew G. Knepley } 743494e7359SMatthew G. Knepley 744f5f57ec0SBarry Smith /*@ 745494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 746494e7359SMatthew G. Knepley 747494e7359SMatthew G. Knepley Not Collective 748494e7359SMatthew G. Knepley 749494e7359SMatthew G. Knepley Input Arguments: 750494e7359SMatthew G. Knepley + dim - The simplex dimension 751a6b92713SMatthew G. Knepley . Nc - The number of components 752dcce0ee2SMatthew G. Knepley . npoints - The number of points in one dimension 753494e7359SMatthew G. Knepley . a - left end of interval (often-1) 754494e7359SMatthew G. Knepley - b - right end of interval (often +1) 755494e7359SMatthew G. Knepley 756744bafbcSMatthew G. Knepley Output Argument: 757552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 758494e7359SMatthew G. Knepley 759494e7359SMatthew G. Knepley Level: intermediate 760494e7359SMatthew G. Knepley 761494e7359SMatthew G. Knepley References: 76296a0c994SBarry Smith . 1. - Karniadakis and Sherwin. FIAT 763494e7359SMatthew G. Knepley 764744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 765494e7359SMatthew G. Knepley @*/ 766dcce0ee2SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 767494e7359SMatthew G. Knepley { 768dcce0ee2SMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints; 769494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 770a6b92713SMatthew G. Knepley PetscInt i, j, k, c; 771494e7359SMatthew G. Knepley PetscErrorCode ierr; 772494e7359SMatthew G. Knepley 773494e7359SMatthew G. Knepley PetscFunctionBegin; 774494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 775dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr); 776dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr); 777494e7359SMatthew G. Knepley switch (dim) { 778707aa5c5SMatthew G. Knepley case 0: 779707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 780707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 781785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 782a6b92713SMatthew G. Knepley ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 783707aa5c5SMatthew G. Knepley x[0] = 0.0; 784a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[c] = 1.0; 785707aa5c5SMatthew G. Knepley break; 786494e7359SMatthew G. Knepley case 1: 787dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr); 788dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr); 789dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i]; 790a6b92713SMatthew G. Knepley ierr = PetscFree(wx);CHKERRQ(ierr); 791494e7359SMatthew G. Knepley break; 792494e7359SMatthew G. Knepley case 2: 793dcce0ee2SMatthew G. Knepley ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr); 794dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 795dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 796dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 797dcce0ee2SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 798dcce0ee2SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 799dcce0ee2SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j]; 800494e7359SMatthew G. Knepley } 801494e7359SMatthew G. Knepley } 802494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 803494e7359SMatthew G. Knepley break; 804494e7359SMatthew G. Knepley case 3: 805dcce0ee2SMatthew G. Knepley ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr); 806dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 807dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 808dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 809dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 810dcce0ee2SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 811dcce0ee2SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 812dcce0ee2SMatthew G. Knepley ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*npoints+j)*npoints+k)*3+0], &x[((i*npoints+j)*npoints+k)*3+1], &x[((i*npoints+j)*npoints+k)*3+2]);CHKERRQ(ierr); 813dcce0ee2SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k]; 814494e7359SMatthew G. Knepley } 815494e7359SMatthew G. Knepley } 816494e7359SMatthew G. Knepley } 817494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 818494e7359SMatthew G. Knepley break; 819494e7359SMatthew G. Knepley default: 820494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 821494e7359SMatthew G. Knepley } 82221454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 8232f5fb066SToby Isaac ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 824dcce0ee2SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 825494e7359SMatthew G. Knepley PetscFunctionReturn(0); 826494e7359SMatthew G. Knepley } 827494e7359SMatthew G. Knepley 828f5f57ec0SBarry Smith /*@ 829b3c0f97bSTom Klotz PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 830b3c0f97bSTom Klotz 831b3c0f97bSTom Klotz Not Collective 832b3c0f97bSTom Klotz 833b3c0f97bSTom Klotz Input Arguments: 834b3c0f97bSTom Klotz + dim - The cell dimension 835b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l 836b3c0f97bSTom Klotz . a - left end of interval (often-1) 837b3c0f97bSTom Klotz - b - right end of interval (often +1) 838b3c0f97bSTom Klotz 839b3c0f97bSTom Klotz Output Argument: 840b3c0f97bSTom Klotz . q - A PetscQuadrature object 841b3c0f97bSTom Klotz 842b3c0f97bSTom Klotz Level: intermediate 843b3c0f97bSTom Klotz 844b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature() 845b3c0f97bSTom Klotz @*/ 846b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 847b3c0f97bSTom Klotz { 848b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 849b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 850b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 851b3c0f97bSTom Klotz const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 852d84b4d08SMatthew G. Knepley PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 853b3c0f97bSTom Klotz PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 854b3c0f97bSTom Klotz PetscReal *x, *w; 855b3c0f97bSTom Klotz PetscInt K, k, npoints; 856b3c0f97bSTom Klotz PetscErrorCode ierr; 857b3c0f97bSTom Klotz 858b3c0f97bSTom Klotz PetscFunctionBegin; 859b3c0f97bSTom Klotz if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 860b3c0f97bSTom Klotz if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 861b3c0f97bSTom Klotz /* Find K such that the weights are < 32 digits of precision */ 862b3c0f97bSTom Klotz for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 8639add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 864b3c0f97bSTom Klotz } 865b3c0f97bSTom Klotz ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 866b3c0f97bSTom Klotz ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 867b3c0f97bSTom Klotz npoints = 2*K-1; 868b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 869b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 870b3c0f97bSTom Klotz /* Center term */ 871b3c0f97bSTom Klotz x[0] = beta; 872b3c0f97bSTom Klotz w[0] = 0.5*alpha*PETSC_PI; 873b3c0f97bSTom Klotz for (k = 1; k < K; ++k) { 8749add2064SThomas Klotz wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 8751118d4bcSLisandro Dalcin xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h)); 876b3c0f97bSTom Klotz x[2*k-1] = -alpha*xk+beta; 877b3c0f97bSTom Klotz w[2*k-1] = wk; 878b3c0f97bSTom Klotz x[2*k+0] = alpha*xk+beta; 879b3c0f97bSTom Klotz w[2*k+0] = wk; 880b3c0f97bSTom Klotz } 881a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr); 882b3c0f97bSTom Klotz PetscFunctionReturn(0); 883b3c0f97bSTom Klotz } 884b3c0f97bSTom Klotz 885b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 886b3c0f97bSTom Klotz { 887b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 888b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 889b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 890b3c0f97bSTom Klotz PetscReal h = 1.0; /* Step size, length between x_k */ 891b3c0f97bSTom Klotz PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 892b3c0f97bSTom Klotz PetscReal osum = 0.0; /* Integral on last level */ 893b3c0f97bSTom Klotz PetscReal psum = 0.0; /* Integral on the level before the last level */ 894b3c0f97bSTom Klotz PetscReal sum; /* Integral on current level */ 895446c295cSMatthew G. Knepley PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 896b3c0f97bSTom Klotz PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 897b3c0f97bSTom Klotz PetscReal wk; /* Quadrature weight at x_k */ 898b3c0f97bSTom Klotz PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 899b3c0f97bSTom Klotz PetscInt d; /* Digits of precision in the integral */ 900b3c0f97bSTom Klotz 901b3c0f97bSTom Klotz PetscFunctionBegin; 902b3c0f97bSTom Klotz if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 903b3c0f97bSTom Klotz /* Center term */ 904b3c0f97bSTom Klotz func(beta, &lval); 905b3c0f97bSTom Klotz sum = 0.5*alpha*PETSC_PI*lval; 906b3c0f97bSTom Klotz /* */ 907b3c0f97bSTom Klotz do { 908b3c0f97bSTom Klotz PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 909b3c0f97bSTom Klotz PetscInt k = 1; 910b3c0f97bSTom Klotz 911b3c0f97bSTom Klotz ++l; 912b3c0f97bSTom Klotz /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 913b3c0f97bSTom Klotz /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 914b3c0f97bSTom Klotz psum = osum; 915b3c0f97bSTom Klotz osum = sum; 916b3c0f97bSTom Klotz h *= 0.5; 917b3c0f97bSTom Klotz sum *= 0.5; 918b3c0f97bSTom Klotz do { 9199add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 920446c295cSMatthew G. Knepley yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 921446c295cSMatthew G. Knepley lx = -alpha*(1.0 - yk)+beta; 922446c295cSMatthew G. Knepley rx = alpha*(1.0 - yk)+beta; 923b3c0f97bSTom Klotz func(lx, &lval); 924b3c0f97bSTom Klotz func(rx, &rval); 925b3c0f97bSTom Klotz lterm = alpha*wk*lval; 926b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 927b3c0f97bSTom Klotz sum += lterm; 928b3c0f97bSTom Klotz rterm = alpha*wk*rval; 929b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 930b3c0f97bSTom Klotz sum += rterm; 931b3c0f97bSTom Klotz ++k; 932b3c0f97bSTom Klotz /* Only need to evaluate every other point on refined levels */ 933b3c0f97bSTom Klotz if (l != 1) ++k; 9349add2064SThomas Klotz } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 935b3c0f97bSTom Klotz 936b3c0f97bSTom Klotz d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 937b3c0f97bSTom Klotz d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 938b3c0f97bSTom Klotz d3 = PetscLog10Real(maxTerm) - p; 93909d48545SBarry Smith if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 94009d48545SBarry Smith else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 941b3c0f97bSTom Klotz d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 9429add2064SThomas Klotz } while (d < digits && l < 12); 943b3c0f97bSTom Klotz *sol = sum; 944e510cb1fSThomas Klotz 945b3c0f97bSTom Klotz PetscFunctionReturn(0); 946b3c0f97bSTom Klotz } 947b3c0f97bSTom Klotz 948*497880caSRichard Tran Mills #if defined(PETSC_HAVE_MPFR) 94929f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 95029f144ccSMatthew G. Knepley { 951e510cb1fSThomas Klotz const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 95229f144ccSMatthew G. Knepley PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 95329f144ccSMatthew G. Knepley mpfr_t alpha; /* Half-width of the integration interval */ 95429f144ccSMatthew G. Knepley mpfr_t beta; /* Center of the integration interval */ 95529f144ccSMatthew G. Knepley mpfr_t h; /* Step size, length between x_k */ 95629f144ccSMatthew G. Knepley mpfr_t osum; /* Integral on last level */ 95729f144ccSMatthew G. Knepley mpfr_t psum; /* Integral on the level before the last level */ 95829f144ccSMatthew G. Knepley mpfr_t sum; /* Integral on current level */ 95929f144ccSMatthew G. Knepley mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 96029f144ccSMatthew G. Knepley mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 96129f144ccSMatthew G. Knepley mpfr_t wk; /* Quadrature weight at x_k */ 96229f144ccSMatthew G. Knepley PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 96329f144ccSMatthew G. Knepley PetscInt d; /* Digits of precision in the integral */ 96429f144ccSMatthew G. Knepley mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 96529f144ccSMatthew G. Knepley 96629f144ccSMatthew G. Knepley PetscFunctionBegin; 96729f144ccSMatthew G. Knepley if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 96829f144ccSMatthew G. Knepley /* Create high precision storage */ 969c9f744b5SMatthew 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); 97029f144ccSMatthew G. Knepley /* Initialization */ 97129f144ccSMatthew G. Knepley mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 97229f144ccSMatthew G. Knepley mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 97329f144ccSMatthew G. Knepley mpfr_set_d(osum, 0.0, MPFR_RNDN); 97429f144ccSMatthew G. Knepley mpfr_set_d(psum, 0.0, MPFR_RNDN); 97529f144ccSMatthew G. Knepley mpfr_set_d(h, 1.0, MPFR_RNDN); 97629f144ccSMatthew G. Knepley mpfr_const_pi(pi2, MPFR_RNDN); 97729f144ccSMatthew G. Knepley mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 97829f144ccSMatthew G. Knepley /* Center term */ 97929f144ccSMatthew G. Knepley func(0.5*(b+a), &lval); 98029f144ccSMatthew G. Knepley mpfr_set(sum, pi2, MPFR_RNDN); 98129f144ccSMatthew G. Knepley mpfr_mul(sum, sum, alpha, MPFR_RNDN); 98229f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 98329f144ccSMatthew G. Knepley /* */ 98429f144ccSMatthew G. Knepley do { 98529f144ccSMatthew G. Knepley PetscReal d1, d2, d3, d4; 98629f144ccSMatthew G. Knepley PetscInt k = 1; 98729f144ccSMatthew G. Knepley 98829f144ccSMatthew G. Knepley ++l; 98929f144ccSMatthew G. Knepley mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 99029f144ccSMatthew G. Knepley /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 99129f144ccSMatthew G. Knepley /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 99229f144ccSMatthew G. Knepley mpfr_set(psum, osum, MPFR_RNDN); 99329f144ccSMatthew G. Knepley mpfr_set(osum, sum, MPFR_RNDN); 99429f144ccSMatthew G. Knepley mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 99529f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 99629f144ccSMatthew G. Knepley do { 99729f144ccSMatthew G. Knepley mpfr_set_si(kh, k, MPFR_RNDN); 99829f144ccSMatthew G. Knepley mpfr_mul(kh, kh, h, MPFR_RNDN); 99929f144ccSMatthew G. Knepley /* Weight */ 100029f144ccSMatthew G. Knepley mpfr_set(wk, h, MPFR_RNDN); 100129f144ccSMatthew G. Knepley mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 100229f144ccSMatthew G. Knepley mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 100329f144ccSMatthew G. Knepley mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 100429f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 100529f144ccSMatthew G. Knepley mpfr_sqr(tmp, tmp, MPFR_RNDN); 100629f144ccSMatthew G. Knepley mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 100729f144ccSMatthew G. Knepley mpfr_div(wk, wk, tmp, MPFR_RNDN); 100829f144ccSMatthew G. Knepley /* Abscissa */ 100929f144ccSMatthew G. Knepley mpfr_set_d(yk, 1.0, MPFR_RNDZ); 101029f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 101129f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 101229f144ccSMatthew G. Knepley mpfr_exp(tmp, msinh, MPFR_RNDN); 101329f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 101429f144ccSMatthew G. Knepley /* Quadrature points */ 101529f144ccSMatthew G. Knepley mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 101629f144ccSMatthew G. Knepley mpfr_mul(lx, lx, alpha, MPFR_RNDU); 101729f144ccSMatthew G. Knepley mpfr_add(lx, lx, beta, MPFR_RNDU); 101829f144ccSMatthew G. Knepley mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 101929f144ccSMatthew G. Knepley mpfr_mul(rx, rx, alpha, MPFR_RNDD); 102029f144ccSMatthew G. Knepley mpfr_add(rx, rx, beta, MPFR_RNDD); 102129f144ccSMatthew G. Knepley /* Evaluation */ 102229f144ccSMatthew G. Knepley func(mpfr_get_d(lx, MPFR_RNDU), &lval); 102329f144ccSMatthew G. Knepley func(mpfr_get_d(rx, MPFR_RNDD), &rval); 102429f144ccSMatthew G. Knepley /* Update */ 102529f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 102629f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 102729f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 102829f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 102929f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 103029f144ccSMatthew G. Knepley mpfr_set(curTerm, tmp, MPFR_RNDN); 103129f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 103229f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 103329f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 103429f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 103529f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 103629f144ccSMatthew G. Knepley mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 103729f144ccSMatthew G. Knepley ++k; 103829f144ccSMatthew G. Knepley /* Only need to evaluate every other point on refined levels */ 103929f144ccSMatthew G. Knepley if (l != 1) ++k; 104029f144ccSMatthew G. Knepley mpfr_log10(tmp, wk, MPFR_RNDN); 104129f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 1042c9f744b5SMatthew G. Knepley } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 104329f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, osum, MPFR_RNDN); 104429f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 104529f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 104629f144ccSMatthew G. Knepley d1 = mpfr_get_d(tmp, MPFR_RNDN); 104729f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, psum, MPFR_RNDN); 104829f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 104929f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 105029f144ccSMatthew G. Knepley d2 = mpfr_get_d(tmp, MPFR_RNDN); 105129f144ccSMatthew G. Knepley mpfr_log10(tmp, maxTerm, MPFR_RNDN); 1052c9f744b5SMatthew G. Knepley d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 105329f144ccSMatthew G. Knepley mpfr_log10(tmp, curTerm, MPFR_RNDN); 105429f144ccSMatthew G. Knepley d4 = mpfr_get_d(tmp, MPFR_RNDN); 105529f144ccSMatthew G. Knepley d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1056b0649871SThomas Klotz } while (d < digits && l < 8); 105729f144ccSMatthew G. Knepley *sol = mpfr_get_d(sum, MPFR_RNDN); 105829f144ccSMatthew G. Knepley /* Cleanup */ 105929f144ccSMatthew G. Knepley mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 106029f144ccSMatthew G. Knepley PetscFunctionReturn(0); 106129f144ccSMatthew G. Knepley } 1062d525116cSMatthew G. Knepley #else 1063fbfcfee5SBarry Smith 1064d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1065d525116cSMatthew G. Knepley { 1066d525116cSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 1067d525116cSMatthew G. Knepley } 106829f144ccSMatthew G. Knepley #endif 106929f144ccSMatthew G. Knepley 1070194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 1071194825f6SJed Brown * A in column-major format 1072194825f6SJed Brown * Ainv in row-major format 1073194825f6SJed Brown * tau has length m 1074194825f6SJed Brown * worksize must be >= max(1,n) 1075194825f6SJed Brown */ 1076194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1077194825f6SJed Brown { 1078194825f6SJed Brown PetscErrorCode ierr; 1079194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1080194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 1081194825f6SJed Brown 1082194825f6SJed Brown PetscFunctionBegin; 1083194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1084194825f6SJed Brown { 1085194825f6SJed Brown PetscInt i,j; 1086dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 1087194825f6SJed Brown for (j=0; j<n; j++) { 1088194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 1089194825f6SJed Brown } 1090194825f6SJed Brown mstride = m; 1091194825f6SJed Brown } 1092194825f6SJed Brown #else 1093194825f6SJed Brown A = A_in; 1094194825f6SJed Brown Ainv = Ainv_out; 1095194825f6SJed Brown #endif 1096194825f6SJed Brown 1097194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 1098194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 1099194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 1100194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 1101194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1102001a771dSBarry Smith PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1103194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 1104194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1105194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1106194825f6SJed Brown 1107194825f6SJed Brown /* Extract an explicit representation of Q */ 1108194825f6SJed Brown Q = Ainv; 1109194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 1110194825f6SJed Brown K = N; /* full rank */ 1111c964aadfSJose E. Roman PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1112194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1113194825f6SJed Brown 1114194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1115194825f6SJed Brown Alpha = 1.0; 1116194825f6SJed Brown ldb = lda; 1117001a771dSBarry Smith PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 1118194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 1119194825f6SJed Brown 1120194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1121194825f6SJed Brown { 1122194825f6SJed Brown PetscInt i; 1123194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 1124194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 1125194825f6SJed Brown } 1126194825f6SJed Brown #endif 1127194825f6SJed Brown PetscFunctionReturn(0); 1128194825f6SJed Brown } 1129194825f6SJed Brown 1130194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 1131194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 1132194825f6SJed Brown { 1133194825f6SJed Brown PetscErrorCode ierr; 1134194825f6SJed Brown PetscReal *Bv; 1135194825f6SJed Brown PetscInt i,j; 1136194825f6SJed Brown 1137194825f6SJed Brown PetscFunctionBegin; 1138785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 1139194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 1140194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 1141194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 1142194825f6SJed Brown for (i=0; i<ninterval; i++) { 1143194825f6SJed Brown for (j=0; j<ndegree; j++) { 1144194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1145194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1146194825f6SJed Brown } 1147194825f6SJed Brown } 1148194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 1149194825f6SJed Brown PetscFunctionReturn(0); 1150194825f6SJed Brown } 1151194825f6SJed Brown 1152194825f6SJed Brown /*@ 1153194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 1154194825f6SJed Brown 1155194825f6SJed Brown Not Collective 1156194825f6SJed Brown 1157194825f6SJed Brown Input Arguments: 1158194825f6SJed Brown + degree - degree of reconstruction polynomial 1159194825f6SJed Brown . nsource - number of source intervals 1160194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 1161194825f6SJed Brown . ntarget - number of target intervals 1162194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 1163194825f6SJed Brown 1164194825f6SJed Brown Output Arguments: 1165194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 1166194825f6SJed Brown 1167194825f6SJed Brown Level: advanced 1168194825f6SJed Brown 1169194825f6SJed Brown .seealso: PetscDTLegendreEval() 1170194825f6SJed Brown @*/ 1171194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 1172194825f6SJed Brown { 1173194825f6SJed Brown PetscErrorCode ierr; 1174194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 1175194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 1176194825f6SJed Brown PetscScalar *tau,*work; 1177194825f6SJed Brown 1178194825f6SJed Brown PetscFunctionBegin; 1179194825f6SJed Brown PetscValidRealPointer(sourcex,3); 1180194825f6SJed Brown PetscValidRealPointer(targetx,5); 1181194825f6SJed Brown PetscValidRealPointer(R,6); 1182194825f6SJed 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); 1183194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 1184194825f6SJed Brown for (i=0; i<nsource; i++) { 118557622a8eSBarry 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]); 1186194825f6SJed Brown } 1187194825f6SJed Brown for (i=0; i<ntarget; i++) { 118857622a8eSBarry 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]); 1189194825f6SJed Brown } 1190194825f6SJed Brown #endif 1191194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 1192194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 1193194825f6SJed Brown center = (xmin + xmax)/2; 1194194825f6SJed Brown hscale = (xmax - xmin)/2; 1195194825f6SJed Brown worksize = nsource; 1196dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 1197dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 1198194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 1199194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 1200194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 1201194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 1202194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 1203194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 1204194825f6SJed Brown for (i=0; i<ntarget; i++) { 1205194825f6SJed Brown PetscReal rowsum = 0; 1206194825f6SJed Brown for (j=0; j<nsource; j++) { 1207194825f6SJed Brown PetscReal sum = 0; 1208194825f6SJed Brown for (k=0; k<degree+1; k++) { 1209194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 1210194825f6SJed Brown } 1211194825f6SJed Brown R[i*nsource+j] = sum; 1212194825f6SJed Brown rowsum += sum; 1213194825f6SJed Brown } 1214194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 1215194825f6SJed Brown } 1216194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 1217194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 1218194825f6SJed Brown PetscFunctionReturn(0); 1219194825f6SJed Brown } 1220