137045ce4SJed Brown /* Discretization tools */ 237045ce4SJed Brown 30c35b76eSJed Brown #include <petscdt.h> /*I "petscdt.h" I*/ 437045ce4SJed Brown #include <petscblaslapack.h> 5af0996ceSBarry Smith #include <petsc/private/petscimpl.h> 6af0996ceSBarry Smith #include <petsc/private/dtimpl.h> 7665c2dedSJed Brown #include <petscviewer.h> 859804f93SMatthew G. Knepley #include <petscdmplex.h> 959804f93SMatthew G. Knepley #include <petscdmshell.h> 1037045ce4SJed Brown 1198c04793SMatthew G. Knepley #if defined(PETSC_HAVE_MPFR) 1298c04793SMatthew G. Knepley #include <mpfr.h> 1398c04793SMatthew G. Knepley #endif 1498c04793SMatthew G. Knepley 150bfcf5a5SMatthew G. Knepley static PetscBool GaussCite = PETSC_FALSE; 160bfcf5a5SMatthew G. Knepley const char GaussCitation[] = "@article{GolubWelsch1969,\n" 170bfcf5a5SMatthew G. Knepley " author = {Golub and Welsch},\n" 180bfcf5a5SMatthew G. Knepley " title = {Calculation of Quadrature Rules},\n" 190bfcf5a5SMatthew G. Knepley " journal = {Math. Comp.},\n" 200bfcf5a5SMatthew G. Knepley " volume = {23},\n" 210bfcf5a5SMatthew G. Knepley " number = {106},\n" 220bfcf5a5SMatthew G. Knepley " pages = {221--230},\n" 230bfcf5a5SMatthew G. Knepley " year = {1969}\n}\n"; 240bfcf5a5SMatthew G. Knepley 25*2cd22861SMatthew G. Knepley 26*2cd22861SMatthew G. Knepley PetscClassId PETSCQUADRATURE_CLASSID = 0; 27*2cd22861SMatthew G. Knepley 2840d8ff71SMatthew G. Knepley /*@ 2940d8ff71SMatthew G. Knepley PetscQuadratureCreate - Create a PetscQuadrature object 3040d8ff71SMatthew G. Knepley 31d083f849SBarry Smith Collective 3240d8ff71SMatthew G. Knepley 3340d8ff71SMatthew G. Knepley Input Parameter: 3440d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object 3540d8ff71SMatthew G. Knepley 3640d8ff71SMatthew G. Knepley Output Parameter: 3740d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 3840d8ff71SMatthew G. Knepley 3940d8ff71SMatthew G. Knepley Level: beginner 4040d8ff71SMatthew G. Knepley 4140d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 4240d8ff71SMatthew G. Knepley @*/ 4321454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 4421454ff5SMatthew G. Knepley { 4521454ff5SMatthew G. Knepley PetscErrorCode ierr; 4621454ff5SMatthew G. Knepley 4721454ff5SMatthew G. Knepley PetscFunctionBegin; 4821454ff5SMatthew G. Knepley PetscValidPointer(q, 2); 49*2cd22861SMatthew G. Knepley ierr = DMInitializePackage();CHKERRQ(ierr); 50*2cd22861SMatthew G. Knepley ierr = PetscHeaderCreate(*q,PETSCQUADRATURE_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 5121454ff5SMatthew G. Knepley (*q)->dim = -1; 52a6b92713SMatthew G. Knepley (*q)->Nc = 1; 53bcede257SMatthew G. Knepley (*q)->order = -1; 5421454ff5SMatthew G. Knepley (*q)->numPoints = 0; 5521454ff5SMatthew G. Knepley (*q)->points = NULL; 5621454ff5SMatthew G. Knepley (*q)->weights = NULL; 5721454ff5SMatthew G. Knepley PetscFunctionReturn(0); 5821454ff5SMatthew G. Knepley } 5921454ff5SMatthew G. Knepley 60c9638911SMatthew G. Knepley /*@ 61c9638911SMatthew G. Knepley PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 62c9638911SMatthew G. Knepley 63d083f849SBarry Smith Collective on q 64c9638911SMatthew G. Knepley 65c9638911SMatthew G. Knepley Input Parameter: 66c9638911SMatthew G. Knepley . q - The PetscQuadrature object 67c9638911SMatthew G. Knepley 68c9638911SMatthew G. Knepley Output Parameter: 69c9638911SMatthew G. Knepley . r - The new PetscQuadrature object 70c9638911SMatthew G. Knepley 71c9638911SMatthew G. Knepley Level: beginner 72c9638911SMatthew G. Knepley 73c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 74c9638911SMatthew G. Knepley @*/ 75c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 76c9638911SMatthew G. Knepley { 77a6b92713SMatthew G. Knepley PetscInt order, dim, Nc, Nq; 78c9638911SMatthew G. Knepley const PetscReal *points, *weights; 79c9638911SMatthew G. Knepley PetscReal *p, *w; 80c9638911SMatthew G. Knepley PetscErrorCode ierr; 81c9638911SMatthew G. Knepley 82c9638911SMatthew G. Knepley PetscFunctionBegin; 83c9638911SMatthew G. Knepley PetscValidPointer(q, 2); 84c9638911SMatthew G. Knepley ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 85c9638911SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 86c9638911SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 87a6b92713SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 88c9638911SMatthew G. Knepley ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 89f0a0bfafSMatthew G. Knepley ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr); 90580bdb30SBarry Smith ierr = PetscArraycpy(p, points, Nq*dim);CHKERRQ(ierr); 91580bdb30SBarry Smith ierr = PetscArraycpy(w, weights, Nc * Nq);CHKERRQ(ierr); 92a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr); 93c9638911SMatthew G. Knepley PetscFunctionReturn(0); 94c9638911SMatthew G. Knepley } 95c9638911SMatthew G. Knepley 9640d8ff71SMatthew G. Knepley /*@ 9740d8ff71SMatthew G. Knepley PetscQuadratureDestroy - Destroys a PetscQuadrature object 9840d8ff71SMatthew G. Knepley 99d083f849SBarry Smith Collective on q 10040d8ff71SMatthew G. Knepley 10140d8ff71SMatthew G. Knepley Input Parameter: 10240d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 10340d8ff71SMatthew G. Knepley 10440d8ff71SMatthew G. Knepley Level: beginner 10540d8ff71SMatthew G. Knepley 10640d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 10740d8ff71SMatthew G. Knepley @*/ 108bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 109bfa639d9SMatthew G. Knepley { 110bfa639d9SMatthew G. Knepley PetscErrorCode ierr; 111bfa639d9SMatthew G. Knepley 112bfa639d9SMatthew G. Knepley PetscFunctionBegin; 11321454ff5SMatthew G. Knepley if (!*q) PetscFunctionReturn(0); 114*2cd22861SMatthew G. Knepley PetscValidHeaderSpecific((*q),PETSCQUADRATURE_CLASSID,1); 11521454ff5SMatthew G. Knepley if (--((PetscObject)(*q))->refct > 0) { 11621454ff5SMatthew G. Knepley *q = NULL; 11721454ff5SMatthew G. Knepley PetscFunctionReturn(0); 11821454ff5SMatthew G. Knepley } 11921454ff5SMatthew G. Knepley ierr = PetscFree((*q)->points);CHKERRQ(ierr); 12021454ff5SMatthew G. Knepley ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 12121454ff5SMatthew G. Knepley ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 12221454ff5SMatthew G. Knepley PetscFunctionReturn(0); 12321454ff5SMatthew G. Knepley } 12421454ff5SMatthew G. Knepley 125bcede257SMatthew G. Knepley /*@ 126a6b92713SMatthew G. Knepley PetscQuadratureGetOrder - Return the order of the method 127bcede257SMatthew G. Knepley 128bcede257SMatthew G. Knepley Not collective 129bcede257SMatthew G. Knepley 130bcede257SMatthew G. Knepley Input Parameter: 131bcede257SMatthew G. Knepley . q - The PetscQuadrature object 132bcede257SMatthew G. Knepley 133bcede257SMatthew G. Knepley Output Parameter: 134bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 135bcede257SMatthew G. Knepley 136bcede257SMatthew G. Knepley Level: intermediate 137bcede257SMatthew G. Knepley 138bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 139bcede257SMatthew G. Knepley @*/ 140bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 141bcede257SMatthew G. Knepley { 142bcede257SMatthew G. Knepley PetscFunctionBegin; 143*2cd22861SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 144bcede257SMatthew G. Knepley PetscValidPointer(order, 2); 145bcede257SMatthew G. Knepley *order = q->order; 146bcede257SMatthew G. Knepley PetscFunctionReturn(0); 147bcede257SMatthew G. Knepley } 148bcede257SMatthew G. Knepley 149bcede257SMatthew G. Knepley /*@ 150a6b92713SMatthew G. Knepley PetscQuadratureSetOrder - Return the order of the method 151bcede257SMatthew G. Knepley 152bcede257SMatthew G. Knepley Not collective 153bcede257SMatthew G. Knepley 154bcede257SMatthew G. Knepley Input Parameters: 155bcede257SMatthew G. Knepley + q - The PetscQuadrature object 156bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 157bcede257SMatthew G. Knepley 158bcede257SMatthew G. Knepley Level: intermediate 159bcede257SMatthew G. Knepley 160bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 161bcede257SMatthew G. Knepley @*/ 162bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 163bcede257SMatthew G. Knepley { 164bcede257SMatthew G. Knepley PetscFunctionBegin; 165*2cd22861SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 166bcede257SMatthew G. Knepley q->order = order; 167bcede257SMatthew G. Knepley PetscFunctionReturn(0); 168bcede257SMatthew G. Knepley } 169bcede257SMatthew G. Knepley 170a6b92713SMatthew G. Knepley /*@ 171a6b92713SMatthew G. Knepley PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated 172a6b92713SMatthew G. Knepley 173a6b92713SMatthew G. Knepley Not collective 174a6b92713SMatthew G. Knepley 175a6b92713SMatthew G. Knepley Input Parameter: 176a6b92713SMatthew G. Knepley . q - The PetscQuadrature object 177a6b92713SMatthew G. Knepley 178a6b92713SMatthew G. Knepley Output Parameter: 179a6b92713SMatthew G. Knepley . Nc - The number of components 180a6b92713SMatthew G. Knepley 181a6b92713SMatthew G. Knepley Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 182a6b92713SMatthew G. Knepley 183a6b92713SMatthew G. Knepley Level: intermediate 184a6b92713SMatthew G. Knepley 185a6b92713SMatthew G. Knepley .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 186a6b92713SMatthew G. Knepley @*/ 187a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc) 188a6b92713SMatthew G. Knepley { 189a6b92713SMatthew G. Knepley PetscFunctionBegin; 190*2cd22861SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 191a6b92713SMatthew G. Knepley PetscValidPointer(Nc, 2); 192a6b92713SMatthew G. Knepley *Nc = q->Nc; 193a6b92713SMatthew G. Knepley PetscFunctionReturn(0); 194a6b92713SMatthew G. Knepley } 195a6b92713SMatthew G. Knepley 196a6b92713SMatthew G. Knepley /*@ 197a6b92713SMatthew G. Knepley PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated 198a6b92713SMatthew G. Knepley 199a6b92713SMatthew G. Knepley Not collective 200a6b92713SMatthew G. Knepley 201a6b92713SMatthew G. Knepley Input Parameters: 202a6b92713SMatthew G. Knepley + q - The PetscQuadrature object 203a6b92713SMatthew G. Knepley - Nc - The number of components 204a6b92713SMatthew G. Knepley 205a6b92713SMatthew G. Knepley Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 206a6b92713SMatthew G. Knepley 207a6b92713SMatthew G. Knepley Level: intermediate 208a6b92713SMatthew G. Knepley 209a6b92713SMatthew G. Knepley .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 210a6b92713SMatthew G. Knepley @*/ 211a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc) 212a6b92713SMatthew G. Knepley { 213a6b92713SMatthew G. Knepley PetscFunctionBegin; 214*2cd22861SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 215a6b92713SMatthew G. Knepley q->Nc = Nc; 216a6b92713SMatthew G. Knepley PetscFunctionReturn(0); 217a6b92713SMatthew G. Knepley } 218a6b92713SMatthew G. Knepley 21940d8ff71SMatthew G. Knepley /*@C 22040d8ff71SMatthew G. Knepley PetscQuadratureGetData - Returns the data defining the quadrature 22140d8ff71SMatthew G. Knepley 22240d8ff71SMatthew G. Knepley Not collective 22340d8ff71SMatthew G. Knepley 22440d8ff71SMatthew G. Knepley Input Parameter: 22540d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 22640d8ff71SMatthew G. Knepley 22740d8ff71SMatthew G. Knepley Output Parameters: 22840d8ff71SMatthew G. Knepley + dim - The spatial dimension 229805e7170SToby Isaac . Nc - The number of components 23040d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 23140d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 23240d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 23340d8ff71SMatthew G. Knepley 23440d8ff71SMatthew G. Knepley Level: intermediate 23540d8ff71SMatthew G. Knepley 23695452b02SPatrick Sanan Fortran Notes: 23795452b02SPatrick Sanan From Fortran you must call PetscQuadratureRestoreData() when you are done with the data 2381fd49c25SBarry Smith 23940d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 24040d8ff71SMatthew G. Knepley @*/ 241a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 24221454ff5SMatthew G. Knepley { 24321454ff5SMatthew G. Knepley PetscFunctionBegin; 244*2cd22861SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 24521454ff5SMatthew G. Knepley if (dim) { 24621454ff5SMatthew G. Knepley PetscValidPointer(dim, 2); 24721454ff5SMatthew G. Knepley *dim = q->dim; 24821454ff5SMatthew G. Knepley } 249a6b92713SMatthew G. Knepley if (Nc) { 250a6b92713SMatthew G. Knepley PetscValidPointer(Nc, 3); 251a6b92713SMatthew G. Knepley *Nc = q->Nc; 252a6b92713SMatthew G. Knepley } 25321454ff5SMatthew G. Knepley if (npoints) { 254a6b92713SMatthew G. Knepley PetscValidPointer(npoints, 4); 25521454ff5SMatthew G. Knepley *npoints = q->numPoints; 25621454ff5SMatthew G. Knepley } 25721454ff5SMatthew G. Knepley if (points) { 258a6b92713SMatthew G. Knepley PetscValidPointer(points, 5); 25921454ff5SMatthew G. Knepley *points = q->points; 26021454ff5SMatthew G. Knepley } 26121454ff5SMatthew G. Knepley if (weights) { 262a6b92713SMatthew G. Knepley PetscValidPointer(weights, 6); 26321454ff5SMatthew G. Knepley *weights = q->weights; 26421454ff5SMatthew G. Knepley } 26521454ff5SMatthew G. Knepley PetscFunctionReturn(0); 26621454ff5SMatthew G. Knepley } 26721454ff5SMatthew G. Knepley 26840d8ff71SMatthew G. Knepley /*@C 26940d8ff71SMatthew G. Knepley PetscQuadratureSetData - Sets the data defining the quadrature 27040d8ff71SMatthew G. Knepley 27140d8ff71SMatthew G. Knepley Not collective 27240d8ff71SMatthew G. Knepley 27340d8ff71SMatthew G. Knepley Input Parameters: 27440d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 27540d8ff71SMatthew G. Knepley . dim - The spatial dimension 276e2b35d93SBarry Smith . Nc - The number of components 27740d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 27840d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 27940d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 28040d8ff71SMatthew G. Knepley 281c99e0549SMatthew G. Knepley Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them. 282f2fd9e53SMatthew G. Knepley 28340d8ff71SMatthew G. Knepley Level: intermediate 28440d8ff71SMatthew G. Knepley 28540d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 28640d8ff71SMatthew G. Knepley @*/ 287a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 28821454ff5SMatthew G. Knepley { 28921454ff5SMatthew G. Knepley PetscFunctionBegin; 290*2cd22861SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 29121454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 292a6b92713SMatthew G. Knepley if (Nc >= 0) q->Nc = Nc; 29321454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 29421454ff5SMatthew G. Knepley if (points) { 29521454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 29621454ff5SMatthew G. Knepley q->points = points; 29721454ff5SMatthew G. Knepley } 29821454ff5SMatthew G. Knepley if (weights) { 29921454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 30021454ff5SMatthew G. Knepley q->weights = weights; 30121454ff5SMatthew G. Knepley } 302f9fd7fdbSMatthew G. Knepley PetscFunctionReturn(0); 303f9fd7fdbSMatthew G. Knepley } 304f9fd7fdbSMatthew G. Knepley 305d9bac1caSLisandro Dalcin static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v) 306d9bac1caSLisandro Dalcin { 307d9bac1caSLisandro Dalcin PetscInt q, d, c; 308d9bac1caSLisandro Dalcin PetscViewerFormat format; 309d9bac1caSLisandro Dalcin PetscErrorCode ierr; 310d9bac1caSLisandro Dalcin 311d9bac1caSLisandro Dalcin PetscFunctionBegin; 312c74b4a09SMatthew G. Knepley if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D) with %D components\n", quad->order, quad->numPoints, quad->dim, quad->Nc);CHKERRQ(ierr);} 313c74b4a09SMatthew G. Knepley else {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);} 314d9bac1caSLisandro Dalcin ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 315d9bac1caSLisandro Dalcin if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0); 316d9bac1caSLisandro Dalcin for (q = 0; q < quad->numPoints; ++q) { 317c74b4a09SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr); 318d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr); 319d9bac1caSLisandro Dalcin for (d = 0; d < quad->dim; ++d) { 320d9bac1caSLisandro Dalcin if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 321d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 322d9bac1caSLisandro Dalcin } 323d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr); 324c74b4a09SMatthew G. Knepley if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);} 325d9bac1caSLisandro Dalcin for (c = 0; c < quad->Nc; ++c) { 326d9bac1caSLisandro Dalcin if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 327c74b4a09SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr); 328d9bac1caSLisandro Dalcin } 329d9bac1caSLisandro Dalcin if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);} 330d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr); 331d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr); 332d9bac1caSLisandro Dalcin } 333d9bac1caSLisandro Dalcin PetscFunctionReturn(0); 334d9bac1caSLisandro Dalcin } 335d9bac1caSLisandro Dalcin 33640d8ff71SMatthew G. Knepley /*@C 33740d8ff71SMatthew G. Knepley PetscQuadratureView - Views a PetscQuadrature object 33840d8ff71SMatthew G. Knepley 339d083f849SBarry Smith Collective on quad 34040d8ff71SMatthew G. Knepley 34140d8ff71SMatthew G. Knepley Input Parameters: 342d9bac1caSLisandro Dalcin + quad - The PetscQuadrature object 34340d8ff71SMatthew G. Knepley - viewer - The PetscViewer object 34440d8ff71SMatthew G. Knepley 34540d8ff71SMatthew G. Knepley Level: beginner 34640d8ff71SMatthew G. Knepley 34740d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 34840d8ff71SMatthew G. Knepley @*/ 349f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 350f9fd7fdbSMatthew G. Knepley { 351d9bac1caSLisandro Dalcin PetscBool iascii; 352f9fd7fdbSMatthew G. Knepley PetscErrorCode ierr; 353f9fd7fdbSMatthew G. Knepley 354f9fd7fdbSMatthew G. Knepley PetscFunctionBegin; 355d9bac1caSLisandro Dalcin PetscValidHeader(quad, 1); 356d9bac1caSLisandro Dalcin if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 357d9bac1caSLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);} 358d9bac1caSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 359d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 360d9bac1caSLisandro Dalcin if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);} 361d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 362bfa639d9SMatthew G. Knepley PetscFunctionReturn(0); 363bfa639d9SMatthew G. Knepley } 364bfa639d9SMatthew G. Knepley 36589710940SMatthew G. Knepley /*@C 36689710940SMatthew G. Knepley PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 36789710940SMatthew G. Knepley 36889710940SMatthew G. Knepley Not collective 36989710940SMatthew G. Knepley 37089710940SMatthew G. Knepley Input Parameter: 37189710940SMatthew G. Knepley + q - The original PetscQuadrature 37289710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into 37389710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement 37489710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement 37589710940SMatthew G. Knepley 37689710940SMatthew G. Knepley Output Parameters: 37789710940SMatthew G. Knepley . dim - The dimension 37889710940SMatthew G. Knepley 37989710940SMatthew G. Knepley Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 38089710940SMatthew G. Knepley 381f5f57ec0SBarry Smith Not available from Fortran 382f5f57ec0SBarry Smith 38389710940SMatthew G. Knepley Level: intermediate 38489710940SMatthew G. Knepley 38589710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 38689710940SMatthew G. Knepley @*/ 38789710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 38889710940SMatthew G. Knepley { 38989710940SMatthew G. Knepley const PetscReal *points, *weights; 39089710940SMatthew G. Knepley PetscReal *pointsRef, *weightsRef; 391a6b92713SMatthew G. Knepley PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e; 39289710940SMatthew G. Knepley PetscErrorCode ierr; 39389710940SMatthew G. Knepley 39489710940SMatthew G. Knepley PetscFunctionBegin; 395*2cd22861SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 39689710940SMatthew G. Knepley PetscValidPointer(v0, 3); 39789710940SMatthew G. Knepley PetscValidPointer(jac, 4); 39889710940SMatthew G. Knepley PetscValidPointer(qref, 5); 39989710940SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 40089710940SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 401a6b92713SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 40289710940SMatthew G. Knepley npointsRef = npoints*numSubelements; 40389710940SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 404a6b92713SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr); 40589710940SMatthew G. Knepley for (c = 0; c < numSubelements; ++c) { 40689710940SMatthew G. Knepley for (p = 0; p < npoints; ++p) { 40789710940SMatthew G. Knepley for (d = 0; d < dim; ++d) { 40889710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 40989710940SMatthew G. Knepley for (e = 0; e < dim; ++e) { 41089710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 41189710940SMatthew G. Knepley } 41289710940SMatthew G. Knepley } 41389710940SMatthew G. Knepley /* Could also use detJ here */ 414a6b92713SMatthew G. Knepley for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements; 41589710940SMatthew G. Knepley } 41689710940SMatthew G. Knepley } 41789710940SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 418a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 41989710940SMatthew G. Knepley PetscFunctionReturn(0); 42089710940SMatthew G. Knepley } 42189710940SMatthew G. Knepley 42237045ce4SJed Brown /*@ 42337045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 42437045ce4SJed Brown 42537045ce4SJed Brown Not Collective 42637045ce4SJed Brown 42737045ce4SJed Brown Input Arguments: 42837045ce4SJed Brown + npoints - number of spatial points to evaluate at 42937045ce4SJed Brown . points - array of locations to evaluate at 43037045ce4SJed Brown . ndegree - number of basis degrees to evaluate 43137045ce4SJed Brown - degrees - sorted array of degrees to evaluate 43237045ce4SJed Brown 43337045ce4SJed Brown Output Arguments: 4340298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 4350298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 4360298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 43737045ce4SJed Brown 43837045ce4SJed Brown Level: intermediate 43937045ce4SJed Brown 44037045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 44137045ce4SJed Brown @*/ 44237045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 44337045ce4SJed Brown { 44437045ce4SJed Brown PetscInt i,maxdegree; 44537045ce4SJed Brown 44637045ce4SJed Brown PetscFunctionBegin; 44737045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 44837045ce4SJed Brown maxdegree = degrees[ndegree-1]; 44937045ce4SJed Brown for (i=0; i<npoints; i++) { 45037045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 45137045ce4SJed Brown PetscInt j,k; 45237045ce4SJed Brown x = points[i]; 45337045ce4SJed Brown pm2 = 0; 45437045ce4SJed Brown pm1 = 1; 45537045ce4SJed Brown pd2 = 0; 45637045ce4SJed Brown pd1 = 0; 45737045ce4SJed Brown pdd2 = 0; 45837045ce4SJed Brown pdd1 = 0; 45937045ce4SJed Brown k = 0; 46037045ce4SJed Brown if (degrees[k] == 0) { 46137045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 46237045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 46337045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 46437045ce4SJed Brown k++; 46537045ce4SJed Brown } 46637045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 46737045ce4SJed Brown PetscReal p,d,dd; 46837045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 46937045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 47037045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 47137045ce4SJed Brown pm2 = pm1; 47237045ce4SJed Brown pm1 = p; 47337045ce4SJed Brown pd2 = pd1; 47437045ce4SJed Brown pd1 = d; 47537045ce4SJed Brown pdd2 = pdd1; 47637045ce4SJed Brown pdd1 = dd; 47737045ce4SJed Brown if (degrees[k] == j) { 47837045ce4SJed Brown if (B) B[i*ndegree+k] = p; 47937045ce4SJed Brown if (D) D[i*ndegree+k] = d; 48037045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 48137045ce4SJed Brown } 48237045ce4SJed Brown } 48337045ce4SJed Brown } 48437045ce4SJed Brown PetscFunctionReturn(0); 48537045ce4SJed Brown } 48637045ce4SJed Brown 48737045ce4SJed Brown /*@ 48837045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 48937045ce4SJed Brown 49037045ce4SJed Brown Not Collective 49137045ce4SJed Brown 49237045ce4SJed Brown Input Arguments: 49337045ce4SJed Brown + npoints - number of points 49437045ce4SJed Brown . a - left end of interval (often-1) 49537045ce4SJed Brown - b - right end of interval (often +1) 49637045ce4SJed Brown 49737045ce4SJed Brown Output Arguments: 49837045ce4SJed Brown + x - quadrature points 49937045ce4SJed Brown - w - quadrature weights 50037045ce4SJed Brown 50137045ce4SJed Brown Level: intermediate 50237045ce4SJed Brown 50337045ce4SJed Brown References: 50496a0c994SBarry Smith . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 50537045ce4SJed Brown 50637045ce4SJed Brown .seealso: PetscDTLegendreEval() 50737045ce4SJed Brown @*/ 50837045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 50937045ce4SJed Brown { 51037045ce4SJed Brown PetscErrorCode ierr; 51137045ce4SJed Brown PetscInt i; 51237045ce4SJed Brown PetscReal *work; 51337045ce4SJed Brown PetscScalar *Z; 51437045ce4SJed Brown PetscBLASInt N,LDZ,info; 51537045ce4SJed Brown 51637045ce4SJed Brown PetscFunctionBegin; 5170bfcf5a5SMatthew G. Knepley ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 51837045ce4SJed Brown /* Set up the Golub-Welsch system */ 51937045ce4SJed Brown for (i=0; i<npoints; i++) { 52037045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 52137045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 52237045ce4SJed Brown } 523dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 524c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 52537045ce4SJed Brown LDZ = N; 52637045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 5278b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 52837045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 5291c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 53037045ce4SJed Brown 53137045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 53237045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 53337045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 53419a57d60SBarry Smith if (x[i] == -0.0) x[i] = 0.0; 53537045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 5360d644c17SKarl Rupp 53788393a60SJed Brown w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 53837045ce4SJed Brown } 53937045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 54037045ce4SJed Brown PetscFunctionReturn(0); 54137045ce4SJed Brown } 542194825f6SJed Brown 5438272889dSSatish Balay static void qAndLEvaluation(PetscInt n, PetscReal x, PetscReal *q, PetscReal *qp, PetscReal *Ln) 5448272889dSSatish Balay /* 5458272889dSSatish Balay Compute the polynomial q(x) = L_{N+1}(x) - L_{n-1}(x) and its derivative in 5468272889dSSatish Balay addition to L_N(x) as these are needed for computing the GLL points via Newton's method. 5478272889dSSatish Balay Reference: "Implementing Spectral Methods for Partial Differential Equations: Algorithms 5488272889dSSatish Balay for Scientists and Engineers" by David A. Kopriva. 5498272889dSSatish Balay */ 5508272889dSSatish Balay { 5518272889dSSatish Balay PetscInt k; 5528272889dSSatish Balay 5538272889dSSatish Balay PetscReal Lnp; 5548272889dSSatish Balay PetscReal Lnp1, Lnp1p; 5558272889dSSatish Balay PetscReal Lnm1, Lnm1p; 5568272889dSSatish Balay PetscReal Lnm2, Lnm2p; 5578272889dSSatish Balay 5588272889dSSatish Balay Lnm1 = 1.0; 5598272889dSSatish Balay *Ln = x; 5608272889dSSatish Balay Lnm1p = 0.0; 5618272889dSSatish Balay Lnp = 1.0; 5628272889dSSatish Balay 5638272889dSSatish Balay for (k=2; k<=n; ++k) { 5648272889dSSatish Balay Lnm2 = Lnm1; 5658272889dSSatish Balay Lnm1 = *Ln; 5668272889dSSatish Balay Lnm2p = Lnm1p; 5678272889dSSatish Balay Lnm1p = Lnp; 5688272889dSSatish Balay *Ln = (2.*((PetscReal)k)-1.)/(1.0*((PetscReal)k))*x*Lnm1 - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm2; 5698272889dSSatish Balay Lnp = Lnm2p + (2.0*((PetscReal)k)-1.)*Lnm1; 5708272889dSSatish Balay } 5718272889dSSatish Balay k = n+1; 5728272889dSSatish Balay Lnp1 = (2.*((PetscReal)k)-1.)/(((PetscReal)k))*x*(*Ln) - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm1; 5738272889dSSatish Balay Lnp1p = Lnm1p + (2.0*((PetscReal)k)-1.)*(*Ln); 5748272889dSSatish Balay *q = Lnp1 - Lnm1; 5758272889dSSatish Balay *qp = Lnp1p - Lnm1p; 5768272889dSSatish Balay } 5778272889dSSatish Balay 5788272889dSSatish Balay /*@C 5798272889dSSatish Balay PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre 5808272889dSSatish Balay nodes of a given size on the domain [-1,1] 5818272889dSSatish Balay 5828272889dSSatish Balay Not Collective 5838272889dSSatish Balay 5848272889dSSatish Balay Input Parameter: 5858272889dSSatish Balay + n - number of grid nodes 586f2e8fe4dShannah_mairs - type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 5878272889dSSatish Balay 5888272889dSSatish Balay Output Arguments: 5898272889dSSatish Balay + x - quadrature points 5908272889dSSatish Balay - w - quadrature weights 5918272889dSSatish Balay 5928272889dSSatish Balay Notes: 5938272889dSSatish Balay For n > 30 the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not 5948272889dSSatish Balay close enough to the desired solution 5958272889dSSatish Balay 5968272889dSSatish Balay These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes 5978272889dSSatish Balay 598a8d69d7bSBarry Smith See https://epubs.siam.org/doi/abs/10.1137/110855442 https://epubs.siam.org/doi/abs/10.1137/120889873 for better ways to compute GLL nodes 5998272889dSSatish Balay 6008272889dSSatish Balay Level: intermediate 6018272889dSSatish Balay 6028272889dSSatish Balay .seealso: PetscDTGaussQuadrature() 6038272889dSSatish Balay 6048272889dSSatish Balay @*/ 605916e780bShannah_mairs PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w) 6068272889dSSatish Balay { 6078272889dSSatish Balay PetscErrorCode ierr; 6088272889dSSatish Balay 6098272889dSSatish Balay PetscFunctionBegin; 6108272889dSSatish Balay if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element"); 6118272889dSSatish Balay 612f2e8fe4dShannah_mairs if (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA) { 6138272889dSSatish Balay PetscReal *M,si; 6148272889dSSatish Balay PetscBLASInt bn,lierr; 6158272889dSSatish Balay PetscReal x0,z0,z1,z2; 6168272889dSSatish Balay PetscInt i,p = npoints - 1,nn; 6178272889dSSatish Balay 6188272889dSSatish Balay x[0] =-1.0; 6198272889dSSatish Balay x[npoints-1] = 1.0; 6208272889dSSatish Balay if (npoints-2 > 0){ 6218272889dSSatish Balay ierr = PetscMalloc1(npoints-1,&M);CHKERRQ(ierr); 6228272889dSSatish Balay for (i=0; i<npoints-2; i++) { 6238272889dSSatish Balay si = ((PetscReal)i)+1.0; 6248272889dSSatish Balay M[i]=0.5*PetscSqrtReal(si*(si+2.0)/((si+0.5)*(si+1.5))); 6258272889dSSatish Balay } 6268272889dSSatish Balay ierr = PetscBLASIntCast(npoints-2,&bn);CHKERRQ(ierr); 627580bdb30SBarry Smith ierr = PetscArrayzero(&x[1],bn);CHKERRQ(ierr); 6288272889dSSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 6298272889dSSatish Balay x0=0; 6308272889dSSatish Balay PetscStackCallBLAS("LAPACKsteqr",LAPACKREALsteqr_("N",&bn,&x[1],M,&x0,&bn,M,&lierr)); 6318272889dSSatish Balay if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in STERF Lapack routine %d",(int)lierr); 6328272889dSSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 6338272889dSSatish Balay ierr = PetscFree(M);CHKERRQ(ierr); 6348272889dSSatish Balay } 6358272889dSSatish Balay if ((npoints-1)%2==0) { 6368272889dSSatish Balay x[(npoints-1)/2] = 0.0; /* hard wire to exactly 0.0 since linear algebra produces nonzero */ 6378272889dSSatish Balay } 6388272889dSSatish Balay 6398272889dSSatish Balay w[0] = w[p] = 2.0/(((PetscReal)(p))*(((PetscReal)p)+1.0)); 6408272889dSSatish Balay z2 = -1.; /* Dummy value to avoid -Wmaybe-initialized */ 6418272889dSSatish Balay for (i=1; i<p; i++) { 6428272889dSSatish Balay x0 = x[i]; 6438272889dSSatish Balay z0 = 1.0; 6448272889dSSatish Balay z1 = x0; 6458272889dSSatish Balay for (nn=1; nn<p; nn++) { 6468272889dSSatish Balay z2 = x0*z1*(2.0*((PetscReal)nn)+1.0)/(((PetscReal)nn)+1.0)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.0)); 6478272889dSSatish Balay z0 = z1; 6488272889dSSatish Balay z1 = z2; 6498272889dSSatish Balay } 6508272889dSSatish Balay w[i]=2.0/(((PetscReal)p)*(((PetscReal)p)+1.0)*z2*z2); 6518272889dSSatish Balay } 6528272889dSSatish Balay } else { 6538272889dSSatish Balay PetscInt j,m; 6548272889dSSatish Balay PetscReal z1,z,q,qp,Ln; 6558272889dSSatish Balay PetscReal *pt; 6568272889dSSatish Balay ierr = PetscMalloc1(npoints,&pt);CHKERRQ(ierr); 6578272889dSSatish Balay 658d410ae54Shannah_mairs if (npoints > 30) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON produces incorrect answers for n > 30"); 6598272889dSSatish Balay x[0] = -1.0; 6608272889dSSatish Balay x[npoints-1] = 1.0; 661feb237baSPierre Jolivet w[0] = w[npoints-1] = 2./(((PetscReal)npoints)*(((PetscReal)npoints)-1.0)); 6628272889dSSatish Balay m = (npoints-1)/2; /* The roots are symmetric, so we only find half of them. */ 6638272889dSSatish Balay for (j=1; j<=m; j++) { /* Loop over the desired roots. */ 6648272889dSSatish Balay z = -1.0*PetscCosReal((PETSC_PI*((PetscReal)j)+0.25)/(((PetscReal)npoints)-1.0))-(3.0/(8.0*(((PetscReal)npoints)-1.0)*PETSC_PI))*(1.0/(((PetscReal)j)+0.25)); 6658272889dSSatish Balay /* Starting with the above approximation to the ith root, we enter */ 6668272889dSSatish Balay /* the main loop of refinement by Newton's method. */ 6678272889dSSatish Balay do { 6688272889dSSatish Balay qAndLEvaluation(npoints-1,z,&q,&qp,&Ln); 6698272889dSSatish Balay z1 = z; 6708272889dSSatish Balay z = z1-q/qp; /* Newton's method. */ 6718272889dSSatish Balay } while (PetscAbs(z-z1) > 10.*PETSC_MACHINE_EPSILON); 6728272889dSSatish Balay qAndLEvaluation(npoints-1,z,&q,&qp,&Ln); 6738272889dSSatish Balay 6748272889dSSatish Balay x[j] = z; 6758272889dSSatish Balay x[npoints-1-j] = -z; /* and put in its symmetric counterpart. */ 6768272889dSSatish Balay w[j] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln); /* Compute the weight */ 6778272889dSSatish Balay w[npoints-1-j] = w[j]; /* and its symmetric counterpart. */ 6788272889dSSatish Balay pt[j]=qp; 6798272889dSSatish Balay } 6808272889dSSatish Balay 6818272889dSSatish Balay if ((npoints-1)%2==0) { 6828272889dSSatish Balay qAndLEvaluation(npoints-1,0.0,&q,&qp,&Ln); 6838272889dSSatish Balay x[(npoints-1)/2] = 0.0; 6848272889dSSatish Balay w[(npoints-1)/2] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln); 6858272889dSSatish Balay } 6868272889dSSatish Balay ierr = PetscFree(pt);CHKERRQ(ierr); 6878272889dSSatish Balay } 6888272889dSSatish Balay PetscFunctionReturn(0); 6898272889dSSatish Balay } 6908272889dSSatish Balay 691744bafbcSMatthew G. Knepley /*@ 692744bafbcSMatthew G. Knepley PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 693744bafbcSMatthew G. Knepley 694744bafbcSMatthew G. Knepley Not Collective 695744bafbcSMatthew G. Knepley 696744bafbcSMatthew G. Knepley Input Arguments: 697744bafbcSMatthew G. Knepley + dim - The spatial dimension 698a6b92713SMatthew G. Knepley . Nc - The number of components 699744bafbcSMatthew G. Knepley . npoints - number of points in one dimension 700744bafbcSMatthew G. Knepley . a - left end of interval (often-1) 701744bafbcSMatthew G. Knepley - b - right end of interval (often +1) 702744bafbcSMatthew G. Knepley 703744bafbcSMatthew G. Knepley Output Argument: 704744bafbcSMatthew G. Knepley . q - A PetscQuadrature object 705744bafbcSMatthew G. Knepley 706744bafbcSMatthew G. Knepley Level: intermediate 707744bafbcSMatthew G. Knepley 708744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 709744bafbcSMatthew G. Knepley @*/ 710a6b92713SMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 711744bafbcSMatthew G. Knepley { 712a6b92713SMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c; 713744bafbcSMatthew G. Knepley PetscReal *x, *w, *xw, *ww; 714744bafbcSMatthew G. Knepley PetscErrorCode ierr; 715744bafbcSMatthew G. Knepley 716744bafbcSMatthew G. Knepley PetscFunctionBegin; 717744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 718a6b92713SMatthew G. Knepley ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr); 719744bafbcSMatthew G. Knepley /* Set up the Golub-Welsch system */ 720744bafbcSMatthew G. Knepley switch (dim) { 721744bafbcSMatthew G. Knepley case 0: 722744bafbcSMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 723744bafbcSMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 724744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 725a6b92713SMatthew G. Knepley ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 726744bafbcSMatthew G. Knepley x[0] = 0.0; 727a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[c] = 1.0; 728744bafbcSMatthew G. Knepley break; 729744bafbcSMatthew G. Knepley case 1: 730a6b92713SMatthew G. Knepley ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr); 731a6b92713SMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr); 732a6b92713SMatthew G. Knepley for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i]; 733a6b92713SMatthew G. Knepley ierr = PetscFree(ww);CHKERRQ(ierr); 734744bafbcSMatthew G. Knepley break; 735744bafbcSMatthew G. Knepley case 2: 736744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 737744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 738744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 739744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 740744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+0] = xw[i]; 741744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+1] = xw[j]; 742a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j]; 743744bafbcSMatthew G. Knepley } 744744bafbcSMatthew G. Knepley } 745744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 746744bafbcSMatthew G. Knepley break; 747744bafbcSMatthew G. Knepley case 3: 748744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 749744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 750744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 751744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 752744bafbcSMatthew G. Knepley for (k = 0; k < npoints; ++k) { 753744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 754744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 755744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 756a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k]; 757744bafbcSMatthew G. Knepley } 758744bafbcSMatthew G. Knepley } 759744bafbcSMatthew G. Knepley } 760744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 761744bafbcSMatthew G. Knepley break; 762744bafbcSMatthew G. Knepley default: 763744bafbcSMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 764744bafbcSMatthew G. Knepley } 765744bafbcSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 7662f5fb066SToby Isaac ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 767a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 768d9bac1caSLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr); 769744bafbcSMatthew G. Knepley PetscFunctionReturn(0); 770744bafbcSMatthew G. Knepley } 771744bafbcSMatthew G. Knepley 772494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 773494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 774494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 775494e7359SMatthew G. Knepley { 776494e7359SMatthew G. Knepley PetscReal f = 1.0; 777494e7359SMatthew G. Knepley PetscInt i; 778494e7359SMatthew G. Knepley 779494e7359SMatthew G. Knepley PetscFunctionBegin; 780494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 781494e7359SMatthew G. Knepley *factorial = f; 782494e7359SMatthew G. Knepley PetscFunctionReturn(0); 783494e7359SMatthew G. Knepley } 784494e7359SMatthew G. Knepley 785494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 786494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 787494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 788494e7359SMatthew G. Knepley { 789494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 790494e7359SMatthew G. Knepley PetscInt k; 791494e7359SMatthew G. Knepley 792494e7359SMatthew G. Knepley PetscFunctionBegin; 793494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 794494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 795494e7359SMatthew G. Knepley apb = a + b; 796494e7359SMatthew G. Knepley pn2 = 1.0; 797494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 798494e7359SMatthew G. Knepley *P = 0.0; 799494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 800494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 801494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 802494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 803494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 804494e7359SMatthew G. Knepley 805494e7359SMatthew G. Knepley a2 = a2 / a1; 806494e7359SMatthew G. Knepley a3 = a3 / a1; 807494e7359SMatthew G. Knepley a4 = a4 / a1; 808494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 809494e7359SMatthew G. Knepley pn2 = pn1; 810494e7359SMatthew G. Knepley pn1 = *P; 811494e7359SMatthew G. Knepley } 812494e7359SMatthew G. Knepley PetscFunctionReturn(0); 813494e7359SMatthew G. Knepley } 814494e7359SMatthew G. Knepley 815494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 816494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 817494e7359SMatthew G. Knepley { 818494e7359SMatthew G. Knepley PetscReal nP; 819494e7359SMatthew G. Knepley PetscErrorCode ierr; 820494e7359SMatthew G. Knepley 821494e7359SMatthew G. Knepley PetscFunctionBegin; 822494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 823494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 824494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 825494e7359SMatthew G. Knepley PetscFunctionReturn(0); 826494e7359SMatthew G. Knepley } 827494e7359SMatthew G. Knepley 828494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 829494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 830494e7359SMatthew G. Knepley { 831494e7359SMatthew G. Knepley PetscFunctionBegin; 832494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 833494e7359SMatthew G. Knepley *eta = y; 834494e7359SMatthew G. Knepley PetscFunctionReturn(0); 835494e7359SMatthew G. Knepley } 836494e7359SMatthew G. Knepley 837494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 838494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 839494e7359SMatthew G. Knepley { 840494e7359SMatthew G. Knepley PetscFunctionBegin; 841494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 842494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 843494e7359SMatthew G. Knepley *zeta = z; 844494e7359SMatthew G. Knepley PetscFunctionReturn(0); 845494e7359SMatthew G. Knepley } 846494e7359SMatthew G. Knepley 847494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 848494e7359SMatthew G. Knepley { 849494e7359SMatthew G. Knepley PetscInt maxIter = 100; 850494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 851a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 852494e7359SMatthew G. Knepley PetscInt k; 853494e7359SMatthew G. Knepley PetscErrorCode ierr; 854494e7359SMatthew G. Knepley 855494e7359SMatthew G. Knepley PetscFunctionBegin; 856a8291ba1SSatish Balay 8578b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 858a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 8590646a658SBarry Smith a2 = PetscTGamma(a + npoints + 1); 8600646a658SBarry Smith a3 = PetscTGamma(b + npoints + 1); 8610646a658SBarry Smith a4 = PetscTGamma(a + b + npoints + 1); 862a8291ba1SSatish Balay #else 86329bcbfd0SToby Isaac { 864d24bbb91SToby Isaac PetscInt ia, ib; 86529bcbfd0SToby Isaac 866d24bbb91SToby Isaac ia = (PetscInt) a; 867d24bbb91SToby Isaac ib = (PetscInt) b; 868d24bbb91SToby 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 */ 869d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr); 870d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr); 871d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr); 87229bcbfd0SToby Isaac } else { 873a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 87429bcbfd0SToby Isaac } 87529bcbfd0SToby Isaac } 876a8291ba1SSatish Balay #endif 877a8291ba1SSatish Balay 878494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 879494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 880494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 881494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 882494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 8838b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 884494e7359SMatthew G. Knepley PetscInt j; 885494e7359SMatthew G. Knepley 886494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 887494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 888494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 889494e7359SMatthew G. Knepley PetscInt i; 890494e7359SMatthew G. Knepley 891494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 892494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 893494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 894494e7359SMatthew G. Knepley delta = f / (fp - f * s); 895494e7359SMatthew G. Knepley r = r - delta; 89677b4d14cSPeter Brune if (PetscAbsReal(delta) < eps) break; 897494e7359SMatthew G. Knepley } 898494e7359SMatthew G. Knepley x[k] = r; 899494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 900494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 901494e7359SMatthew G. Knepley } 902494e7359SMatthew G. Knepley PetscFunctionReturn(0); 903494e7359SMatthew G. Knepley } 904494e7359SMatthew G. Knepley 905f5f57ec0SBarry Smith /*@ 906494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 907494e7359SMatthew G. Knepley 908494e7359SMatthew G. Knepley Not Collective 909494e7359SMatthew G. Knepley 910494e7359SMatthew G. Knepley Input Arguments: 911494e7359SMatthew G. Knepley + dim - The simplex dimension 912a6b92713SMatthew G. Knepley . Nc - The number of components 913dcce0ee2SMatthew G. Knepley . npoints - The number of points in one dimension 914494e7359SMatthew G. Knepley . a - left end of interval (often-1) 915494e7359SMatthew G. Knepley - b - right end of interval (often +1) 916494e7359SMatthew G. Knepley 917744bafbcSMatthew G. Knepley Output Argument: 918552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 919494e7359SMatthew G. Knepley 920494e7359SMatthew G. Knepley Level: intermediate 921494e7359SMatthew G. Knepley 922494e7359SMatthew G. Knepley References: 92396a0c994SBarry Smith . 1. - Karniadakis and Sherwin. FIAT 924494e7359SMatthew G. Knepley 925744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 926494e7359SMatthew G. Knepley @*/ 927dcce0ee2SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 928494e7359SMatthew G. Knepley { 929dcce0ee2SMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints; 930494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 931a6b92713SMatthew G. Knepley PetscInt i, j, k, c; 932494e7359SMatthew G. Knepley PetscErrorCode ierr; 933494e7359SMatthew G. Knepley 934494e7359SMatthew G. Knepley PetscFunctionBegin; 935494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 936dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr); 937dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr); 938494e7359SMatthew G. Knepley switch (dim) { 939707aa5c5SMatthew G. Knepley case 0: 940707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 941707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 942785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 943a6b92713SMatthew G. Knepley ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 944707aa5c5SMatthew G. Knepley x[0] = 0.0; 945a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[c] = 1.0; 946707aa5c5SMatthew G. Knepley break; 947494e7359SMatthew G. Knepley case 1: 948dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr); 949dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr); 950dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i]; 951a6b92713SMatthew G. Knepley ierr = PetscFree(wx);CHKERRQ(ierr); 952494e7359SMatthew G. Knepley break; 953494e7359SMatthew G. Knepley case 2: 954dcce0ee2SMatthew G. Knepley ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr); 955dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 956dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 957dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 958dcce0ee2SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 959dcce0ee2SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 960dcce0ee2SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j]; 961494e7359SMatthew G. Knepley } 962494e7359SMatthew G. Knepley } 963494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 964494e7359SMatthew G. Knepley break; 965494e7359SMatthew G. Knepley case 3: 966dcce0ee2SMatthew G. Knepley ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr); 967dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 968dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 969dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 970dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 971dcce0ee2SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 972dcce0ee2SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 973dcce0ee2SMatthew 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); 974dcce0ee2SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k]; 975494e7359SMatthew G. Knepley } 976494e7359SMatthew G. Knepley } 977494e7359SMatthew G. Knepley } 978494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 979494e7359SMatthew G. Knepley break; 980494e7359SMatthew G. Knepley default: 981494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 982494e7359SMatthew G. Knepley } 98321454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 9842f5fb066SToby Isaac ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 985dcce0ee2SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 986d9bac1caSLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr); 987494e7359SMatthew G. Knepley PetscFunctionReturn(0); 988494e7359SMatthew G. Knepley } 989494e7359SMatthew G. Knepley 990f5f57ec0SBarry Smith /*@ 991b3c0f97bSTom Klotz PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 992b3c0f97bSTom Klotz 993b3c0f97bSTom Klotz Not Collective 994b3c0f97bSTom Klotz 995b3c0f97bSTom Klotz Input Arguments: 996b3c0f97bSTom Klotz + dim - The cell dimension 997b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l 998b3c0f97bSTom Klotz . a - left end of interval (often-1) 999b3c0f97bSTom Klotz - b - right end of interval (often +1) 1000b3c0f97bSTom Klotz 1001b3c0f97bSTom Klotz Output Argument: 1002b3c0f97bSTom Klotz . q - A PetscQuadrature object 1003b3c0f97bSTom Klotz 1004b3c0f97bSTom Klotz Level: intermediate 1005b3c0f97bSTom Klotz 1006b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature() 1007b3c0f97bSTom Klotz @*/ 1008b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 1009b3c0f97bSTom Klotz { 1010b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 1011b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1012b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1013b3c0f97bSTom Klotz const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 1014d84b4d08SMatthew G. Knepley PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 1015b3c0f97bSTom Klotz PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 1016b3c0f97bSTom Klotz PetscReal *x, *w; 1017b3c0f97bSTom Klotz PetscInt K, k, npoints; 1018b3c0f97bSTom Klotz PetscErrorCode ierr; 1019b3c0f97bSTom Klotz 1020b3c0f97bSTom Klotz PetscFunctionBegin; 1021b3c0f97bSTom Klotz if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 1022b3c0f97bSTom Klotz if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 1023b3c0f97bSTom Klotz /* Find K such that the weights are < 32 digits of precision */ 1024b3c0f97bSTom Klotz for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 10259add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 1026b3c0f97bSTom Klotz } 1027b3c0f97bSTom Klotz ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1028b3c0f97bSTom Klotz ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 1029b3c0f97bSTom Klotz npoints = 2*K-1; 1030b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 1031b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 1032b3c0f97bSTom Klotz /* Center term */ 1033b3c0f97bSTom Klotz x[0] = beta; 1034b3c0f97bSTom Klotz w[0] = 0.5*alpha*PETSC_PI; 1035b3c0f97bSTom Klotz for (k = 1; k < K; ++k) { 10369add2064SThomas Klotz wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 10371118d4bcSLisandro Dalcin xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h)); 1038b3c0f97bSTom Klotz x[2*k-1] = -alpha*xk+beta; 1039b3c0f97bSTom Klotz w[2*k-1] = wk; 1040b3c0f97bSTom Klotz x[2*k+0] = alpha*xk+beta; 1041b3c0f97bSTom Klotz w[2*k+0] = wk; 1042b3c0f97bSTom Klotz } 1043a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr); 1044b3c0f97bSTom Klotz PetscFunctionReturn(0); 1045b3c0f97bSTom Klotz } 1046b3c0f97bSTom Klotz 1047b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1048b3c0f97bSTom Klotz { 1049b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 1050b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1051b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1052b3c0f97bSTom Klotz PetscReal h = 1.0; /* Step size, length between x_k */ 1053b3c0f97bSTom Klotz PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 1054b3c0f97bSTom Klotz PetscReal osum = 0.0; /* Integral on last level */ 1055b3c0f97bSTom Klotz PetscReal psum = 0.0; /* Integral on the level before the last level */ 1056b3c0f97bSTom Klotz PetscReal sum; /* Integral on current level */ 1057446c295cSMatthew G. Knepley PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 1058b3c0f97bSTom Klotz PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 1059b3c0f97bSTom Klotz PetscReal wk; /* Quadrature weight at x_k */ 1060b3c0f97bSTom Klotz PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 1061b3c0f97bSTom Klotz PetscInt d; /* Digits of precision in the integral */ 1062b3c0f97bSTom Klotz 1063b3c0f97bSTom Klotz PetscFunctionBegin; 1064b3c0f97bSTom Klotz if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 1065b3c0f97bSTom Klotz /* Center term */ 1066b3c0f97bSTom Klotz func(beta, &lval); 1067b3c0f97bSTom Klotz sum = 0.5*alpha*PETSC_PI*lval; 1068b3c0f97bSTom Klotz /* */ 1069b3c0f97bSTom Klotz do { 1070b3c0f97bSTom Klotz PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 1071b3c0f97bSTom Klotz PetscInt k = 1; 1072b3c0f97bSTom Klotz 1073b3c0f97bSTom Klotz ++l; 1074b3c0f97bSTom Klotz /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 1075b3c0f97bSTom Klotz /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 1076b3c0f97bSTom Klotz psum = osum; 1077b3c0f97bSTom Klotz osum = sum; 1078b3c0f97bSTom Klotz h *= 0.5; 1079b3c0f97bSTom Klotz sum *= 0.5; 1080b3c0f97bSTom Klotz do { 10819add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1082446c295cSMatthew G. Knepley yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1083446c295cSMatthew G. Knepley lx = -alpha*(1.0 - yk)+beta; 1084446c295cSMatthew G. Knepley rx = alpha*(1.0 - yk)+beta; 1085b3c0f97bSTom Klotz func(lx, &lval); 1086b3c0f97bSTom Klotz func(rx, &rval); 1087b3c0f97bSTom Klotz lterm = alpha*wk*lval; 1088b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 1089b3c0f97bSTom Klotz sum += lterm; 1090b3c0f97bSTom Klotz rterm = alpha*wk*rval; 1091b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 1092b3c0f97bSTom Klotz sum += rterm; 1093b3c0f97bSTom Klotz ++k; 1094b3c0f97bSTom Klotz /* Only need to evaluate every other point on refined levels */ 1095b3c0f97bSTom Klotz if (l != 1) ++k; 10969add2064SThomas Klotz } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 1097b3c0f97bSTom Klotz 1098b3c0f97bSTom Klotz d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 1099b3c0f97bSTom Klotz d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 1100b3c0f97bSTom Klotz d3 = PetscLog10Real(maxTerm) - p; 110109d48545SBarry Smith if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 110209d48545SBarry Smith else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 1103b3c0f97bSTom Klotz d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 11049add2064SThomas Klotz } while (d < digits && l < 12); 1105b3c0f97bSTom Klotz *sol = sum; 1106e510cb1fSThomas Klotz 1107b3c0f97bSTom Klotz PetscFunctionReturn(0); 1108b3c0f97bSTom Klotz } 1109b3c0f97bSTom Klotz 1110497880caSRichard Tran Mills #if defined(PETSC_HAVE_MPFR) 111129f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 111229f144ccSMatthew G. Knepley { 1113e510cb1fSThomas Klotz const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 111429f144ccSMatthew G. Knepley PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 111529f144ccSMatthew G. Knepley mpfr_t alpha; /* Half-width of the integration interval */ 111629f144ccSMatthew G. Knepley mpfr_t beta; /* Center of the integration interval */ 111729f144ccSMatthew G. Knepley mpfr_t h; /* Step size, length between x_k */ 111829f144ccSMatthew G. Knepley mpfr_t osum; /* Integral on last level */ 111929f144ccSMatthew G. Knepley mpfr_t psum; /* Integral on the level before the last level */ 112029f144ccSMatthew G. Knepley mpfr_t sum; /* Integral on current level */ 112129f144ccSMatthew G. Knepley mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 112229f144ccSMatthew G. Knepley mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 112329f144ccSMatthew G. Knepley mpfr_t wk; /* Quadrature weight at x_k */ 112429f144ccSMatthew G. Knepley PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 112529f144ccSMatthew G. Knepley PetscInt d; /* Digits of precision in the integral */ 112629f144ccSMatthew G. Knepley mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 112729f144ccSMatthew G. Knepley 112829f144ccSMatthew G. Knepley PetscFunctionBegin; 112929f144ccSMatthew G. Knepley if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 113029f144ccSMatthew G. Knepley /* Create high precision storage */ 1131c9f744b5SMatthew 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); 113229f144ccSMatthew G. Knepley /* Initialization */ 113329f144ccSMatthew G. Knepley mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 113429f144ccSMatthew G. Knepley mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 113529f144ccSMatthew G. Knepley mpfr_set_d(osum, 0.0, MPFR_RNDN); 113629f144ccSMatthew G. Knepley mpfr_set_d(psum, 0.0, MPFR_RNDN); 113729f144ccSMatthew G. Knepley mpfr_set_d(h, 1.0, MPFR_RNDN); 113829f144ccSMatthew G. Knepley mpfr_const_pi(pi2, MPFR_RNDN); 113929f144ccSMatthew G. Knepley mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 114029f144ccSMatthew G. Knepley /* Center term */ 114129f144ccSMatthew G. Knepley func(0.5*(b+a), &lval); 114229f144ccSMatthew G. Knepley mpfr_set(sum, pi2, MPFR_RNDN); 114329f144ccSMatthew G. Knepley mpfr_mul(sum, sum, alpha, MPFR_RNDN); 114429f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 114529f144ccSMatthew G. Knepley /* */ 114629f144ccSMatthew G. Knepley do { 114729f144ccSMatthew G. Knepley PetscReal d1, d2, d3, d4; 114829f144ccSMatthew G. Knepley PetscInt k = 1; 114929f144ccSMatthew G. Knepley 115029f144ccSMatthew G. Knepley ++l; 115129f144ccSMatthew G. Knepley mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 115229f144ccSMatthew G. Knepley /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 115329f144ccSMatthew G. Knepley /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 115429f144ccSMatthew G. Knepley mpfr_set(psum, osum, MPFR_RNDN); 115529f144ccSMatthew G. Knepley mpfr_set(osum, sum, MPFR_RNDN); 115629f144ccSMatthew G. Knepley mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 115729f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 115829f144ccSMatthew G. Knepley do { 115929f144ccSMatthew G. Knepley mpfr_set_si(kh, k, MPFR_RNDN); 116029f144ccSMatthew G. Knepley mpfr_mul(kh, kh, h, MPFR_RNDN); 116129f144ccSMatthew G. Knepley /* Weight */ 116229f144ccSMatthew G. Knepley mpfr_set(wk, h, MPFR_RNDN); 116329f144ccSMatthew G. Knepley mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 116429f144ccSMatthew G. Knepley mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 116529f144ccSMatthew G. Knepley mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 116629f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 116729f144ccSMatthew G. Knepley mpfr_sqr(tmp, tmp, MPFR_RNDN); 116829f144ccSMatthew G. Knepley mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 116929f144ccSMatthew G. Knepley mpfr_div(wk, wk, tmp, MPFR_RNDN); 117029f144ccSMatthew G. Knepley /* Abscissa */ 117129f144ccSMatthew G. Knepley mpfr_set_d(yk, 1.0, MPFR_RNDZ); 117229f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 117329f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 117429f144ccSMatthew G. Knepley mpfr_exp(tmp, msinh, MPFR_RNDN); 117529f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 117629f144ccSMatthew G. Knepley /* Quadrature points */ 117729f144ccSMatthew G. Knepley mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 117829f144ccSMatthew G. Knepley mpfr_mul(lx, lx, alpha, MPFR_RNDU); 117929f144ccSMatthew G. Knepley mpfr_add(lx, lx, beta, MPFR_RNDU); 118029f144ccSMatthew G. Knepley mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 118129f144ccSMatthew G. Knepley mpfr_mul(rx, rx, alpha, MPFR_RNDD); 118229f144ccSMatthew G. Knepley mpfr_add(rx, rx, beta, MPFR_RNDD); 118329f144ccSMatthew G. Knepley /* Evaluation */ 118429f144ccSMatthew G. Knepley func(mpfr_get_d(lx, MPFR_RNDU), &lval); 118529f144ccSMatthew G. Knepley func(mpfr_get_d(rx, MPFR_RNDD), &rval); 118629f144ccSMatthew G. Knepley /* Update */ 118729f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 118829f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 118929f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 119029f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 119129f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 119229f144ccSMatthew G. Knepley mpfr_set(curTerm, tmp, MPFR_RNDN); 119329f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 119429f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 119529f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 119629f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 119729f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 119829f144ccSMatthew G. Knepley mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 119929f144ccSMatthew G. Knepley ++k; 120029f144ccSMatthew G. Knepley /* Only need to evaluate every other point on refined levels */ 120129f144ccSMatthew G. Knepley if (l != 1) ++k; 120229f144ccSMatthew G. Knepley mpfr_log10(tmp, wk, MPFR_RNDN); 120329f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 1204c9f744b5SMatthew G. Knepley } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 120529f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, osum, MPFR_RNDN); 120629f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 120729f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 120829f144ccSMatthew G. Knepley d1 = mpfr_get_d(tmp, MPFR_RNDN); 120929f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, psum, MPFR_RNDN); 121029f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 121129f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 121229f144ccSMatthew G. Knepley d2 = mpfr_get_d(tmp, MPFR_RNDN); 121329f144ccSMatthew G. Knepley mpfr_log10(tmp, maxTerm, MPFR_RNDN); 1214c9f744b5SMatthew G. Knepley d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 121529f144ccSMatthew G. Knepley mpfr_log10(tmp, curTerm, MPFR_RNDN); 121629f144ccSMatthew G. Knepley d4 = mpfr_get_d(tmp, MPFR_RNDN); 121729f144ccSMatthew G. Knepley d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1218b0649871SThomas Klotz } while (d < digits && l < 8); 121929f144ccSMatthew G. Knepley *sol = mpfr_get_d(sum, MPFR_RNDN); 122029f144ccSMatthew G. Knepley /* Cleanup */ 122129f144ccSMatthew G. Knepley mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 122229f144ccSMatthew G. Knepley PetscFunctionReturn(0); 122329f144ccSMatthew G. Knepley } 1224d525116cSMatthew G. Knepley #else 1225fbfcfee5SBarry Smith 1226d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1227d525116cSMatthew G. Knepley { 1228d525116cSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 1229d525116cSMatthew G. Knepley } 123029f144ccSMatthew G. Knepley #endif 123129f144ccSMatthew G. Knepley 1232194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 1233194825f6SJed Brown * A in column-major format 1234194825f6SJed Brown * Ainv in row-major format 1235194825f6SJed Brown * tau has length m 1236194825f6SJed Brown * worksize must be >= max(1,n) 1237194825f6SJed Brown */ 1238194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1239194825f6SJed Brown { 1240194825f6SJed Brown PetscErrorCode ierr; 1241194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1242194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 1243194825f6SJed Brown 1244194825f6SJed Brown PetscFunctionBegin; 1245194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1246194825f6SJed Brown { 1247194825f6SJed Brown PetscInt i,j; 1248dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 1249194825f6SJed Brown for (j=0; j<n; j++) { 1250194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 1251194825f6SJed Brown } 1252194825f6SJed Brown mstride = m; 1253194825f6SJed Brown } 1254194825f6SJed Brown #else 1255194825f6SJed Brown A = A_in; 1256194825f6SJed Brown Ainv = Ainv_out; 1257194825f6SJed Brown #endif 1258194825f6SJed Brown 1259194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 1260194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 1261194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 1262194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 1263194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1264001a771dSBarry Smith PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1265194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 1266194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1267194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1268194825f6SJed Brown 1269194825f6SJed Brown /* Extract an explicit representation of Q */ 1270194825f6SJed Brown Q = Ainv; 1271580bdb30SBarry Smith ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr); 1272194825f6SJed Brown K = N; /* full rank */ 1273c964aadfSJose E. Roman PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1274194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1275194825f6SJed Brown 1276194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1277194825f6SJed Brown Alpha = 1.0; 1278194825f6SJed Brown ldb = lda; 1279001a771dSBarry Smith PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 1280194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 1281194825f6SJed Brown 1282194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1283194825f6SJed Brown { 1284194825f6SJed Brown PetscInt i; 1285194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 1286194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 1287194825f6SJed Brown } 1288194825f6SJed Brown #endif 1289194825f6SJed Brown PetscFunctionReturn(0); 1290194825f6SJed Brown } 1291194825f6SJed Brown 1292194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 1293194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 1294194825f6SJed Brown { 1295194825f6SJed Brown PetscErrorCode ierr; 1296194825f6SJed Brown PetscReal *Bv; 1297194825f6SJed Brown PetscInt i,j; 1298194825f6SJed Brown 1299194825f6SJed Brown PetscFunctionBegin; 1300785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 1301194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 1302194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 1303194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 1304194825f6SJed Brown for (i=0; i<ninterval; i++) { 1305194825f6SJed Brown for (j=0; j<ndegree; j++) { 1306194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1307194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1308194825f6SJed Brown } 1309194825f6SJed Brown } 1310194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 1311194825f6SJed Brown PetscFunctionReturn(0); 1312194825f6SJed Brown } 1313194825f6SJed Brown 1314194825f6SJed Brown /*@ 1315194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 1316194825f6SJed Brown 1317194825f6SJed Brown Not Collective 1318194825f6SJed Brown 1319194825f6SJed Brown Input Arguments: 1320194825f6SJed Brown + degree - degree of reconstruction polynomial 1321194825f6SJed Brown . nsource - number of source intervals 1322194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 1323194825f6SJed Brown . ntarget - number of target intervals 1324194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 1325194825f6SJed Brown 1326194825f6SJed Brown Output Arguments: 1327194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 1328194825f6SJed Brown 1329194825f6SJed Brown Level: advanced 1330194825f6SJed Brown 1331194825f6SJed Brown .seealso: PetscDTLegendreEval() 1332194825f6SJed Brown @*/ 1333194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 1334194825f6SJed Brown { 1335194825f6SJed Brown PetscErrorCode ierr; 1336194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 1337194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 1338194825f6SJed Brown PetscScalar *tau,*work; 1339194825f6SJed Brown 1340194825f6SJed Brown PetscFunctionBegin; 1341194825f6SJed Brown PetscValidRealPointer(sourcex,3); 1342194825f6SJed Brown PetscValidRealPointer(targetx,5); 1343194825f6SJed Brown PetscValidRealPointer(R,6); 1344194825f6SJed 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); 1345194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 1346194825f6SJed Brown for (i=0; i<nsource; i++) { 134757622a8eSBarry 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]); 1348194825f6SJed Brown } 1349194825f6SJed Brown for (i=0; i<ntarget; i++) { 135057622a8eSBarry 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]); 1351194825f6SJed Brown } 1352194825f6SJed Brown #endif 1353194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 1354194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 1355194825f6SJed Brown center = (xmin + xmax)/2; 1356194825f6SJed Brown hscale = (xmax - xmin)/2; 1357194825f6SJed Brown worksize = nsource; 1358dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 1359dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 1360194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 1361194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 1362194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 1363194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 1364194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 1365194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 1366194825f6SJed Brown for (i=0; i<ntarget; i++) { 1367194825f6SJed Brown PetscReal rowsum = 0; 1368194825f6SJed Brown for (j=0; j<nsource; j++) { 1369194825f6SJed Brown PetscReal sum = 0; 1370194825f6SJed Brown for (k=0; k<degree+1; k++) { 1371194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 1372194825f6SJed Brown } 1373194825f6SJed Brown R[i*nsource+j] = sum; 1374194825f6SJed Brown rowsum += sum; 1375194825f6SJed Brown } 1376194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 1377194825f6SJed Brown } 1378194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 1379194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 1380194825f6SJed Brown PetscFunctionReturn(0); 1381194825f6SJed Brown } 1382916e780bShannah_mairs 1383916e780bShannah_mairs /*@C 1384916e780bShannah_mairs PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points 1385916e780bShannah_mairs 1386916e780bShannah_mairs Not Collective 1387916e780bShannah_mairs 1388916e780bShannah_mairs Input Parameter: 1389916e780bShannah_mairs + n - the number of GLL nodes 1390916e780bShannah_mairs . nodes - the GLL nodes 1391916e780bShannah_mairs . weights - the GLL weights 1392916e780bShannah_mairs . f - the function values at the nodes 1393916e780bShannah_mairs 1394916e780bShannah_mairs Output Parameter: 1395916e780bShannah_mairs . in - the value of the integral 1396916e780bShannah_mairs 1397916e780bShannah_mairs Level: beginner 1398916e780bShannah_mairs 1399916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature() 1400916e780bShannah_mairs 1401916e780bShannah_mairs @*/ 1402916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in) 1403916e780bShannah_mairs { 1404916e780bShannah_mairs PetscInt i; 1405916e780bShannah_mairs 1406916e780bShannah_mairs PetscFunctionBegin; 1407916e780bShannah_mairs *in = 0.; 1408916e780bShannah_mairs for (i=0; i<n; i++) { 1409916e780bShannah_mairs *in += f[i]*f[i]*weights[i]; 1410916e780bShannah_mairs } 1411916e780bShannah_mairs PetscFunctionReturn(0); 1412916e780bShannah_mairs } 1413916e780bShannah_mairs 1414916e780bShannah_mairs /*@C 1415916e780bShannah_mairs PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element 1416916e780bShannah_mairs 1417916e780bShannah_mairs Not Collective 1418916e780bShannah_mairs 1419916e780bShannah_mairs Input Parameter: 1420916e780bShannah_mairs + n - the number of GLL nodes 1421916e780bShannah_mairs . nodes - the GLL nodes 1422916e780bShannah_mairs . weights - the GLL weights 1423916e780bShannah_mairs 1424916e780bShannah_mairs Output Parameter: 1425916e780bShannah_mairs . A - the stiffness element 1426916e780bShannah_mairs 1427916e780bShannah_mairs Level: beginner 1428916e780bShannah_mairs 1429916e780bShannah_mairs Notes: 1430916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy() 1431916e780bShannah_mairs 1432916e780bShannah_mairs You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented (the array is symmetric) 1433916e780bShannah_mairs 1434916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1435916e780bShannah_mairs 1436916e780bShannah_mairs @*/ 1437916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1438916e780bShannah_mairs { 1439916e780bShannah_mairs PetscReal **A; 1440916e780bShannah_mairs PetscErrorCode ierr; 1441916e780bShannah_mairs const PetscReal *gllnodes = nodes; 1442916e780bShannah_mairs const PetscInt p = n-1; 1443916e780bShannah_mairs PetscReal z0,z1,z2 = -1,x,Lpj,Lpr; 1444916e780bShannah_mairs PetscInt i,j,nn,r; 1445916e780bShannah_mairs 1446916e780bShannah_mairs PetscFunctionBegin; 1447916e780bShannah_mairs ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1448916e780bShannah_mairs ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1449916e780bShannah_mairs for (i=1; i<n; i++) A[i] = A[i-1]+n; 1450916e780bShannah_mairs 1451916e780bShannah_mairs for (j=1; j<p; j++) { 1452916e780bShannah_mairs x = gllnodes[j]; 1453916e780bShannah_mairs z0 = 1.; 1454916e780bShannah_mairs z1 = x; 1455916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1456916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1457916e780bShannah_mairs z0 = z1; 1458916e780bShannah_mairs z1 = z2; 1459916e780bShannah_mairs } 1460916e780bShannah_mairs Lpj=z2; 1461916e780bShannah_mairs for (r=1; r<p; r++) { 1462916e780bShannah_mairs if (r == j) { 1463916e780bShannah_mairs A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj); 1464916e780bShannah_mairs } else { 1465916e780bShannah_mairs x = gllnodes[r]; 1466916e780bShannah_mairs z0 = 1.; 1467916e780bShannah_mairs z1 = x; 1468916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1469916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1470916e780bShannah_mairs z0 = z1; 1471916e780bShannah_mairs z1 = z2; 1472916e780bShannah_mairs } 1473916e780bShannah_mairs Lpr = z2; 1474916e780bShannah_mairs A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r])); 1475916e780bShannah_mairs } 1476916e780bShannah_mairs } 1477916e780bShannah_mairs } 1478916e780bShannah_mairs for (j=1; j<p+1; j++) { 1479916e780bShannah_mairs x = gllnodes[j]; 1480916e780bShannah_mairs z0 = 1.; 1481916e780bShannah_mairs z1 = x; 1482916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1483916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1484916e780bShannah_mairs z0 = z1; 1485916e780bShannah_mairs z1 = z2; 1486916e780bShannah_mairs } 1487916e780bShannah_mairs Lpj = z2; 1488916e780bShannah_mairs A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j])); 1489916e780bShannah_mairs A[0][j] = A[j][0]; 1490916e780bShannah_mairs } 1491916e780bShannah_mairs for (j=0; j<p; j++) { 1492916e780bShannah_mairs x = gllnodes[j]; 1493916e780bShannah_mairs z0 = 1.; 1494916e780bShannah_mairs z1 = x; 1495916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1496916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1497916e780bShannah_mairs z0 = z1; 1498916e780bShannah_mairs z1 = z2; 1499916e780bShannah_mairs } 1500916e780bShannah_mairs Lpj=z2; 1501916e780bShannah_mairs 1502916e780bShannah_mairs A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j])); 1503916e780bShannah_mairs A[j][p] = A[p][j]; 1504916e780bShannah_mairs } 1505916e780bShannah_mairs A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.; 1506916e780bShannah_mairs A[p][p]=A[0][0]; 1507916e780bShannah_mairs *AA = A; 1508916e780bShannah_mairs PetscFunctionReturn(0); 1509916e780bShannah_mairs } 1510916e780bShannah_mairs 1511916e780bShannah_mairs /*@C 1512916e780bShannah_mairs PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element 1513916e780bShannah_mairs 1514916e780bShannah_mairs Not Collective 1515916e780bShannah_mairs 1516916e780bShannah_mairs Input Parameter: 1517916e780bShannah_mairs + n - the number of GLL nodes 1518916e780bShannah_mairs . nodes - the GLL nodes 1519916e780bShannah_mairs . weights - the GLL weightss 1520916e780bShannah_mairs - A - the stiffness element 1521916e780bShannah_mairs 1522916e780bShannah_mairs Level: beginner 1523916e780bShannah_mairs 1524916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate() 1525916e780bShannah_mairs 1526916e780bShannah_mairs @*/ 1527916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1528916e780bShannah_mairs { 1529916e780bShannah_mairs PetscErrorCode ierr; 1530916e780bShannah_mairs 1531916e780bShannah_mairs PetscFunctionBegin; 1532916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1533916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1534916e780bShannah_mairs *AA = NULL; 1535916e780bShannah_mairs PetscFunctionReturn(0); 1536916e780bShannah_mairs } 1537916e780bShannah_mairs 1538916e780bShannah_mairs /*@C 1539916e780bShannah_mairs PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element 1540916e780bShannah_mairs 1541916e780bShannah_mairs Not Collective 1542916e780bShannah_mairs 1543916e780bShannah_mairs Input Parameter: 1544916e780bShannah_mairs + n - the number of GLL nodes 1545916e780bShannah_mairs . nodes - the GLL nodes 1546916e780bShannah_mairs . weights - the GLL weights 1547916e780bShannah_mairs 1548916e780bShannah_mairs Output Parameter: 1549916e780bShannah_mairs . AA - the stiffness element 1550916e780bShannah_mairs - AAT - the transpose of AA (pass in NULL if you do not need this array) 1551916e780bShannah_mairs 1552916e780bShannah_mairs Level: beginner 1553916e780bShannah_mairs 1554916e780bShannah_mairs Notes: 1555916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementGradientDestroy() 1556916e780bShannah_mairs 1557916e780bShannah_mairs You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1558916e780bShannah_mairs 1559916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1560916e780bShannah_mairs 1561916e780bShannah_mairs @*/ 1562916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1563916e780bShannah_mairs { 1564916e780bShannah_mairs PetscReal **A, **AT = NULL; 1565916e780bShannah_mairs PetscErrorCode ierr; 1566916e780bShannah_mairs const PetscReal *gllnodes = nodes; 1567916e780bShannah_mairs const PetscInt p = n-1; 1568916e780bShannah_mairs PetscReal q,qp,Li, Lj,d0; 1569916e780bShannah_mairs PetscInt i,j; 1570916e780bShannah_mairs 1571916e780bShannah_mairs PetscFunctionBegin; 1572916e780bShannah_mairs ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1573916e780bShannah_mairs ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1574916e780bShannah_mairs for (i=1; i<n; i++) A[i] = A[i-1]+n; 1575916e780bShannah_mairs 1576916e780bShannah_mairs if (AAT) { 1577916e780bShannah_mairs ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr); 1578916e780bShannah_mairs ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr); 1579916e780bShannah_mairs for (i=1; i<n; i++) AT[i] = AT[i-1]+n; 1580916e780bShannah_mairs } 1581916e780bShannah_mairs 1582916e780bShannah_mairs if (n==1) {A[0][0] = 0.;} 1583916e780bShannah_mairs d0 = (PetscReal)p*((PetscReal)p+1.)/4.; 1584916e780bShannah_mairs for (i=0; i<n; i++) { 1585916e780bShannah_mairs for (j=0; j<n; j++) { 1586916e780bShannah_mairs A[i][j] = 0.; 1587fdd31e58Shannah_mairs qAndLEvaluation(p,gllnodes[i],&q,&qp,&Li); 1588fdd31e58Shannah_mairs qAndLEvaluation(p,gllnodes[j],&q,&qp,&Lj); 1589916e780bShannah_mairs if (i!=j) A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j])); 1590916e780bShannah_mairs if ((j==i) && (i==0)) A[i][j] = -d0; 1591916e780bShannah_mairs if (j==i && i==p) A[i][j] = d0; 1592916e780bShannah_mairs if (AT) AT[j][i] = A[i][j]; 1593916e780bShannah_mairs } 1594916e780bShannah_mairs } 1595916e780bShannah_mairs if (AAT) *AAT = AT; 1596916e780bShannah_mairs *AA = A; 1597916e780bShannah_mairs PetscFunctionReturn(0); 1598916e780bShannah_mairs } 1599916e780bShannah_mairs 1600916e780bShannah_mairs /*@C 1601916e780bShannah_mairs PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate() 1602916e780bShannah_mairs 1603916e780bShannah_mairs Not Collective 1604916e780bShannah_mairs 1605916e780bShannah_mairs Input Parameter: 1606916e780bShannah_mairs + n - the number of GLL nodes 1607916e780bShannah_mairs . nodes - the GLL nodes 1608916e780bShannah_mairs . weights - the GLL weights 1609916e780bShannah_mairs . AA - the stiffness element 1610916e780bShannah_mairs - AAT - the transpose of the element 1611916e780bShannah_mairs 1612916e780bShannah_mairs Level: beginner 1613916e780bShannah_mairs 1614916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate() 1615916e780bShannah_mairs 1616916e780bShannah_mairs @*/ 1617916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1618916e780bShannah_mairs { 1619916e780bShannah_mairs PetscErrorCode ierr; 1620916e780bShannah_mairs 1621916e780bShannah_mairs PetscFunctionBegin; 1622916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1623916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1624916e780bShannah_mairs *AA = NULL; 1625916e780bShannah_mairs if (*AAT) { 1626916e780bShannah_mairs ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr); 1627916e780bShannah_mairs ierr = PetscFree(*AAT);CHKERRQ(ierr); 1628916e780bShannah_mairs *AAT = NULL; 1629916e780bShannah_mairs } 1630916e780bShannah_mairs PetscFunctionReturn(0); 1631916e780bShannah_mairs } 1632916e780bShannah_mairs 1633916e780bShannah_mairs /*@C 1634916e780bShannah_mairs PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element 1635916e780bShannah_mairs 1636916e780bShannah_mairs Not Collective 1637916e780bShannah_mairs 1638916e780bShannah_mairs Input Parameter: 1639916e780bShannah_mairs + n - the number of GLL nodes 1640916e780bShannah_mairs . nodes - the GLL nodes 1641916e780bShannah_mairs . weights - the GLL weightss 1642916e780bShannah_mairs 1643916e780bShannah_mairs Output Parameter: 1644916e780bShannah_mairs . AA - the stiffness element 1645916e780bShannah_mairs 1646916e780bShannah_mairs Level: beginner 1647916e780bShannah_mairs 1648916e780bShannah_mairs Notes: 1649916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy() 1650916e780bShannah_mairs 1651916e780bShannah_mairs This is the same as the Gradient operator multiplied by the diagonal mass matrix 1652916e780bShannah_mairs 1653916e780bShannah_mairs You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1654916e780bShannah_mairs 1655916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy() 1656916e780bShannah_mairs 1657916e780bShannah_mairs @*/ 1658916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1659916e780bShannah_mairs { 1660916e780bShannah_mairs PetscReal **D; 1661916e780bShannah_mairs PetscErrorCode ierr; 1662916e780bShannah_mairs const PetscReal *gllweights = weights; 1663916e780bShannah_mairs const PetscInt glln = n; 1664916e780bShannah_mairs PetscInt i,j; 1665916e780bShannah_mairs 1666916e780bShannah_mairs PetscFunctionBegin; 1667916e780bShannah_mairs ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr); 1668916e780bShannah_mairs for (i=0; i<glln; i++){ 1669916e780bShannah_mairs for (j=0; j<glln; j++) { 1670916e780bShannah_mairs D[i][j] = gllweights[i]*D[i][j]; 1671916e780bShannah_mairs } 1672916e780bShannah_mairs } 1673916e780bShannah_mairs *AA = D; 1674916e780bShannah_mairs PetscFunctionReturn(0); 1675916e780bShannah_mairs } 1676916e780bShannah_mairs 1677916e780bShannah_mairs /*@C 1678916e780bShannah_mairs PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element 1679916e780bShannah_mairs 1680916e780bShannah_mairs Not Collective 1681916e780bShannah_mairs 1682916e780bShannah_mairs Input Parameter: 1683916e780bShannah_mairs + n - the number of GLL nodes 1684916e780bShannah_mairs . nodes - the GLL nodes 1685916e780bShannah_mairs . weights - the GLL weights 1686916e780bShannah_mairs - A - advection 1687916e780bShannah_mairs 1688916e780bShannah_mairs Level: beginner 1689916e780bShannah_mairs 1690916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate() 1691916e780bShannah_mairs 1692916e780bShannah_mairs @*/ 1693916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1694916e780bShannah_mairs { 1695916e780bShannah_mairs PetscErrorCode ierr; 1696916e780bShannah_mairs 1697916e780bShannah_mairs PetscFunctionBegin; 1698916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1699916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1700916e780bShannah_mairs *AA = NULL; 1701916e780bShannah_mairs PetscFunctionReturn(0); 1702916e780bShannah_mairs } 1703916e780bShannah_mairs 1704916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1705916e780bShannah_mairs { 1706916e780bShannah_mairs PetscReal **A; 1707916e780bShannah_mairs PetscErrorCode ierr; 1708916e780bShannah_mairs const PetscReal *gllweights = weights; 1709916e780bShannah_mairs const PetscInt glln = n; 1710916e780bShannah_mairs PetscInt i,j; 1711916e780bShannah_mairs 1712916e780bShannah_mairs PetscFunctionBegin; 1713916e780bShannah_mairs ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr); 1714916e780bShannah_mairs ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr); 1715916e780bShannah_mairs for (i=1; i<glln; i++) A[i] = A[i-1]+glln; 1716916e780bShannah_mairs if (glln==1) {A[0][0] = 0.;} 1717916e780bShannah_mairs for (i=0; i<glln; i++) { 1718916e780bShannah_mairs for (j=0; j<glln; j++) { 1719916e780bShannah_mairs A[i][j] = 0.; 1720916e780bShannah_mairs if (j==i) A[i][j] = gllweights[i]; 1721916e780bShannah_mairs } 1722916e780bShannah_mairs } 1723916e780bShannah_mairs *AA = A; 1724916e780bShannah_mairs PetscFunctionReturn(0); 1725916e780bShannah_mairs } 1726916e780bShannah_mairs 1727916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1728916e780bShannah_mairs { 1729916e780bShannah_mairs PetscErrorCode ierr; 1730916e780bShannah_mairs 1731916e780bShannah_mairs PetscFunctionBegin; 1732916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1733916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1734916e780bShannah_mairs *AA = NULL; 1735916e780bShannah_mairs PetscFunctionReturn(0); 1736916e780bShannah_mairs } 1737