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 2540d8ff71SMatthew G. Knepley /*@ 2640d8ff71SMatthew G. Knepley PetscQuadratureCreate - Create a PetscQuadrature object 2740d8ff71SMatthew G. Knepley 28*d083f849SBarry Smith Collective 2940d8ff71SMatthew G. Knepley 3040d8ff71SMatthew G. Knepley Input Parameter: 3140d8ff71SMatthew G. Knepley . comm - The communicator for the PetscQuadrature object 3240d8ff71SMatthew G. Knepley 3340d8ff71SMatthew G. Knepley Output Parameter: 3440d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 3540d8ff71SMatthew G. Knepley 3640d8ff71SMatthew G. Knepley Level: beginner 3740d8ff71SMatthew G. Knepley 3840d8ff71SMatthew G. Knepley .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 3940d8ff71SMatthew G. Knepley @*/ 4021454ff5SMatthew G. Knepley PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 4121454ff5SMatthew G. Knepley { 4221454ff5SMatthew G. Knepley PetscErrorCode ierr; 4321454ff5SMatthew G. Knepley 4421454ff5SMatthew G. Knepley PetscFunctionBegin; 4521454ff5SMatthew G. Knepley PetscValidPointer(q, 2); 46623436dcSMatthew G. Knepley ierr = PetscSysInitializePackage();CHKERRQ(ierr); 4773107ff1SLisandro Dalcin ierr = PetscHeaderCreate(*q,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 4821454ff5SMatthew G. Knepley (*q)->dim = -1; 49a6b92713SMatthew G. Knepley (*q)->Nc = 1; 50bcede257SMatthew G. Knepley (*q)->order = -1; 5121454ff5SMatthew G. Knepley (*q)->numPoints = 0; 5221454ff5SMatthew G. Knepley (*q)->points = NULL; 5321454ff5SMatthew G. Knepley (*q)->weights = NULL; 5421454ff5SMatthew G. Knepley PetscFunctionReturn(0); 5521454ff5SMatthew G. Knepley } 5621454ff5SMatthew G. Knepley 57c9638911SMatthew G. Knepley /*@ 58c9638911SMatthew G. Knepley PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 59c9638911SMatthew G. Knepley 60*d083f849SBarry Smith Collective on q 61c9638911SMatthew G. Knepley 62c9638911SMatthew G. Knepley Input Parameter: 63c9638911SMatthew G. Knepley . q - The PetscQuadrature object 64c9638911SMatthew G. Knepley 65c9638911SMatthew G. Knepley Output Parameter: 66c9638911SMatthew G. Knepley . r - The new PetscQuadrature object 67c9638911SMatthew G. Knepley 68c9638911SMatthew G. Knepley Level: beginner 69c9638911SMatthew G. Knepley 70c9638911SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 71c9638911SMatthew G. Knepley @*/ 72c9638911SMatthew G. Knepley PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 73c9638911SMatthew G. Knepley { 74a6b92713SMatthew G. Knepley PetscInt order, dim, Nc, Nq; 75c9638911SMatthew G. Knepley const PetscReal *points, *weights; 76c9638911SMatthew G. Knepley PetscReal *p, *w; 77c9638911SMatthew G. Knepley PetscErrorCode ierr; 78c9638911SMatthew G. Knepley 79c9638911SMatthew G. Knepley PetscFunctionBegin; 80c9638911SMatthew G. Knepley PetscValidPointer(q, 2); 81c9638911SMatthew G. Knepley ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 82c9638911SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 83c9638911SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 84a6b92713SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 85c9638911SMatthew G. Knepley ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 86f0a0bfafSMatthew G. Knepley ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr); 87c9638911SMatthew G. Knepley ierr = PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));CHKERRQ(ierr); 88a6b92713SMatthew G. Knepley ierr = PetscMemcpy(w, weights, Nc * Nq * sizeof(PetscReal));CHKERRQ(ierr); 89a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr); 90c9638911SMatthew G. Knepley PetscFunctionReturn(0); 91c9638911SMatthew G. Knepley } 92c9638911SMatthew G. Knepley 9340d8ff71SMatthew G. Knepley /*@ 9440d8ff71SMatthew G. Knepley PetscQuadratureDestroy - Destroys a PetscQuadrature object 9540d8ff71SMatthew G. Knepley 96*d083f849SBarry Smith Collective on q 9740d8ff71SMatthew G. Knepley 9840d8ff71SMatthew G. Knepley Input Parameter: 9940d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 10040d8ff71SMatthew G. Knepley 10140d8ff71SMatthew G. Knepley Level: beginner 10240d8ff71SMatthew G. Knepley 10340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 10440d8ff71SMatthew G. Knepley @*/ 105bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 106bfa639d9SMatthew G. Knepley { 107bfa639d9SMatthew G. Knepley PetscErrorCode ierr; 108bfa639d9SMatthew G. Knepley 109bfa639d9SMatthew G. Knepley PetscFunctionBegin; 11021454ff5SMatthew G. Knepley if (!*q) PetscFunctionReturn(0); 11121454ff5SMatthew G. Knepley PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 11221454ff5SMatthew G. Knepley if (--((PetscObject)(*q))->refct > 0) { 11321454ff5SMatthew G. Knepley *q = NULL; 11421454ff5SMatthew G. Knepley PetscFunctionReturn(0); 11521454ff5SMatthew G. Knepley } 11621454ff5SMatthew G. Knepley ierr = PetscFree((*q)->points);CHKERRQ(ierr); 11721454ff5SMatthew G. Knepley ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 11821454ff5SMatthew G. Knepley ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 11921454ff5SMatthew G. Knepley PetscFunctionReturn(0); 12021454ff5SMatthew G. Knepley } 12121454ff5SMatthew G. Knepley 122bcede257SMatthew G. Knepley /*@ 123a6b92713SMatthew G. Knepley PetscQuadratureGetOrder - Return the order of the method 124bcede257SMatthew G. Knepley 125bcede257SMatthew G. Knepley Not collective 126bcede257SMatthew G. Knepley 127bcede257SMatthew G. Knepley Input Parameter: 128bcede257SMatthew G. Knepley . q - The PetscQuadrature object 129bcede257SMatthew G. Knepley 130bcede257SMatthew G. Knepley Output Parameter: 131bcede257SMatthew G. Knepley . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 132bcede257SMatthew G. Knepley 133bcede257SMatthew G. Knepley Level: intermediate 134bcede257SMatthew G. Knepley 135bcede257SMatthew G. Knepley .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 136bcede257SMatthew G. Knepley @*/ 137bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 138bcede257SMatthew G. Knepley { 139bcede257SMatthew G. Knepley PetscFunctionBegin; 140bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 141bcede257SMatthew G. Knepley PetscValidPointer(order, 2); 142bcede257SMatthew G. Knepley *order = q->order; 143bcede257SMatthew G. Knepley PetscFunctionReturn(0); 144bcede257SMatthew G. Knepley } 145bcede257SMatthew G. Knepley 146bcede257SMatthew G. Knepley /*@ 147a6b92713SMatthew G. Knepley PetscQuadratureSetOrder - Return the order of the method 148bcede257SMatthew G. Knepley 149bcede257SMatthew G. Knepley Not collective 150bcede257SMatthew G. Knepley 151bcede257SMatthew G. Knepley Input Parameters: 152bcede257SMatthew G. Knepley + q - The PetscQuadrature object 153bcede257SMatthew G. Knepley - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 154bcede257SMatthew G. Knepley 155bcede257SMatthew G. Knepley Level: intermediate 156bcede257SMatthew G. Knepley 157bcede257SMatthew G. Knepley .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 158bcede257SMatthew G. Knepley @*/ 159bcede257SMatthew G. Knepley PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 160bcede257SMatthew G. Knepley { 161bcede257SMatthew G. Knepley PetscFunctionBegin; 162bcede257SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 163bcede257SMatthew G. Knepley q->order = order; 164bcede257SMatthew G. Knepley PetscFunctionReturn(0); 165bcede257SMatthew G. Knepley } 166bcede257SMatthew G. Knepley 167a6b92713SMatthew G. Knepley /*@ 168a6b92713SMatthew G. Knepley PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated 169a6b92713SMatthew G. Knepley 170a6b92713SMatthew G. Knepley Not collective 171a6b92713SMatthew G. Knepley 172a6b92713SMatthew G. Knepley Input Parameter: 173a6b92713SMatthew G. Knepley . q - The PetscQuadrature object 174a6b92713SMatthew G. Knepley 175a6b92713SMatthew G. Knepley Output Parameter: 176a6b92713SMatthew G. Knepley . Nc - The number of components 177a6b92713SMatthew G. Knepley 178a6b92713SMatthew G. Knepley Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 179a6b92713SMatthew G. Knepley 180a6b92713SMatthew G. Knepley Level: intermediate 181a6b92713SMatthew G. Knepley 182a6b92713SMatthew G. Knepley .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 183a6b92713SMatthew G. Knepley @*/ 184a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc) 185a6b92713SMatthew G. Knepley { 186a6b92713SMatthew G. Knepley PetscFunctionBegin; 187a6b92713SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 188a6b92713SMatthew G. Knepley PetscValidPointer(Nc, 2); 189a6b92713SMatthew G. Knepley *Nc = q->Nc; 190a6b92713SMatthew G. Knepley PetscFunctionReturn(0); 191a6b92713SMatthew G. Knepley } 192a6b92713SMatthew G. Knepley 193a6b92713SMatthew G. Knepley /*@ 194a6b92713SMatthew G. Knepley PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated 195a6b92713SMatthew G. Knepley 196a6b92713SMatthew G. Knepley Not collective 197a6b92713SMatthew G. Knepley 198a6b92713SMatthew G. Knepley Input Parameters: 199a6b92713SMatthew G. Knepley + q - The PetscQuadrature object 200a6b92713SMatthew G. Knepley - Nc - The number of components 201a6b92713SMatthew G. Knepley 202a6b92713SMatthew G. Knepley Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 203a6b92713SMatthew G. Knepley 204a6b92713SMatthew G. Knepley Level: intermediate 205a6b92713SMatthew G. Knepley 206a6b92713SMatthew G. Knepley .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 207a6b92713SMatthew G. Knepley @*/ 208a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc) 209a6b92713SMatthew G. Knepley { 210a6b92713SMatthew G. Knepley PetscFunctionBegin; 211a6b92713SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 212a6b92713SMatthew G. Knepley q->Nc = Nc; 213a6b92713SMatthew G. Knepley PetscFunctionReturn(0); 214a6b92713SMatthew G. Knepley } 215a6b92713SMatthew G. Knepley 21640d8ff71SMatthew G. Knepley /*@C 21740d8ff71SMatthew G. Knepley PetscQuadratureGetData - Returns the data defining the quadrature 21840d8ff71SMatthew G. Knepley 21940d8ff71SMatthew G. Knepley Not collective 22040d8ff71SMatthew G. Knepley 22140d8ff71SMatthew G. Knepley Input Parameter: 22240d8ff71SMatthew G. Knepley . q - The PetscQuadrature object 22340d8ff71SMatthew G. Knepley 22440d8ff71SMatthew G. Knepley Output Parameters: 22540d8ff71SMatthew G. Knepley + dim - The spatial dimension 226805e7170SToby Isaac . Nc - The number of components 22740d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 22840d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 22940d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 23040d8ff71SMatthew G. Knepley 23140d8ff71SMatthew G. Knepley Level: intermediate 23240d8ff71SMatthew G. Knepley 23395452b02SPatrick Sanan Fortran Notes: 23495452b02SPatrick Sanan From Fortran you must call PetscQuadratureRestoreData() when you are done with the data 2351fd49c25SBarry Smith 23640d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 23740d8ff71SMatthew G. Knepley @*/ 238a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 23921454ff5SMatthew G. Knepley { 24021454ff5SMatthew G. Knepley PetscFunctionBegin; 24121454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 24221454ff5SMatthew G. Knepley if (dim) { 24321454ff5SMatthew G. Knepley PetscValidPointer(dim, 2); 24421454ff5SMatthew G. Knepley *dim = q->dim; 24521454ff5SMatthew G. Knepley } 246a6b92713SMatthew G. Knepley if (Nc) { 247a6b92713SMatthew G. Knepley PetscValidPointer(Nc, 3); 248a6b92713SMatthew G. Knepley *Nc = q->Nc; 249a6b92713SMatthew G. Knepley } 25021454ff5SMatthew G. Knepley if (npoints) { 251a6b92713SMatthew G. Knepley PetscValidPointer(npoints, 4); 25221454ff5SMatthew G. Knepley *npoints = q->numPoints; 25321454ff5SMatthew G. Knepley } 25421454ff5SMatthew G. Knepley if (points) { 255a6b92713SMatthew G. Knepley PetscValidPointer(points, 5); 25621454ff5SMatthew G. Knepley *points = q->points; 25721454ff5SMatthew G. Knepley } 25821454ff5SMatthew G. Knepley if (weights) { 259a6b92713SMatthew G. Knepley PetscValidPointer(weights, 6); 26021454ff5SMatthew G. Knepley *weights = q->weights; 26121454ff5SMatthew G. Knepley } 26221454ff5SMatthew G. Knepley PetscFunctionReturn(0); 26321454ff5SMatthew G. Knepley } 26421454ff5SMatthew G. Knepley 26540d8ff71SMatthew G. Knepley /*@C 26640d8ff71SMatthew G. Knepley PetscQuadratureSetData - Sets the data defining the quadrature 26740d8ff71SMatthew G. Knepley 26840d8ff71SMatthew G. Knepley Not collective 26940d8ff71SMatthew G. Knepley 27040d8ff71SMatthew G. Knepley Input Parameters: 27140d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 27240d8ff71SMatthew G. Knepley . dim - The spatial dimension 273e2b35d93SBarry Smith . Nc - The number of components 27440d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 27540d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 27640d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 27740d8ff71SMatthew G. Knepley 278c99e0549SMatthew 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. 279f2fd9e53SMatthew G. Knepley 28040d8ff71SMatthew G. Knepley Level: intermediate 28140d8ff71SMatthew G. Knepley 28240d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 28340d8ff71SMatthew G. Knepley @*/ 284a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 28521454ff5SMatthew G. Knepley { 28621454ff5SMatthew G. Knepley PetscFunctionBegin; 28721454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 28821454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 289a6b92713SMatthew G. Knepley if (Nc >= 0) q->Nc = Nc; 29021454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 29121454ff5SMatthew G. Knepley if (points) { 29221454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 29321454ff5SMatthew G. Knepley q->points = points; 29421454ff5SMatthew G. Knepley } 29521454ff5SMatthew G. Knepley if (weights) { 29621454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 29721454ff5SMatthew G. Knepley q->weights = weights; 29821454ff5SMatthew G. Knepley } 299f9fd7fdbSMatthew G. Knepley PetscFunctionReturn(0); 300f9fd7fdbSMatthew G. Knepley } 301f9fd7fdbSMatthew G. Knepley 302d9bac1caSLisandro Dalcin static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v) 303d9bac1caSLisandro Dalcin { 304d9bac1caSLisandro Dalcin PetscInt q, d, c; 305d9bac1caSLisandro Dalcin PetscViewerFormat format; 306d9bac1caSLisandro Dalcin PetscErrorCode ierr; 307d9bac1caSLisandro Dalcin 308d9bac1caSLisandro Dalcin PetscFunctionBegin; 309c74b4a09SMatthew 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);} 310c74b4a09SMatthew G. Knepley else {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);} 311d9bac1caSLisandro Dalcin ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 312d9bac1caSLisandro Dalcin if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0); 313d9bac1caSLisandro Dalcin for (q = 0; q < quad->numPoints; ++q) { 314c74b4a09SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr); 315d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr); 316d9bac1caSLisandro Dalcin for (d = 0; d < quad->dim; ++d) { 317d9bac1caSLisandro Dalcin if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 318d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 319d9bac1caSLisandro Dalcin } 320d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr); 321c74b4a09SMatthew G. Knepley if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);} 322d9bac1caSLisandro Dalcin for (c = 0; c < quad->Nc; ++c) { 323d9bac1caSLisandro Dalcin if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 324c74b4a09SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr); 325d9bac1caSLisandro Dalcin } 326d9bac1caSLisandro Dalcin if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);} 327d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr); 328d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr); 329d9bac1caSLisandro Dalcin } 330d9bac1caSLisandro Dalcin PetscFunctionReturn(0); 331d9bac1caSLisandro Dalcin } 332d9bac1caSLisandro Dalcin 33340d8ff71SMatthew G. Knepley /*@C 33440d8ff71SMatthew G. Knepley PetscQuadratureView - Views a PetscQuadrature object 33540d8ff71SMatthew G. Knepley 336*d083f849SBarry Smith Collective on quad 33740d8ff71SMatthew G. Knepley 33840d8ff71SMatthew G. Knepley Input Parameters: 339d9bac1caSLisandro Dalcin + quad - The PetscQuadrature object 34040d8ff71SMatthew G. Knepley - viewer - The PetscViewer object 34140d8ff71SMatthew G. Knepley 34240d8ff71SMatthew G. Knepley Level: beginner 34340d8ff71SMatthew G. Knepley 34440d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 34540d8ff71SMatthew G. Knepley @*/ 346f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 347f9fd7fdbSMatthew G. Knepley { 348d9bac1caSLisandro Dalcin PetscBool iascii; 349f9fd7fdbSMatthew G. Knepley PetscErrorCode ierr; 350f9fd7fdbSMatthew G. Knepley 351f9fd7fdbSMatthew G. Knepley PetscFunctionBegin; 352d9bac1caSLisandro Dalcin PetscValidHeader(quad, 1); 353d9bac1caSLisandro Dalcin if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 354d9bac1caSLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);} 355d9bac1caSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 356d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 357d9bac1caSLisandro Dalcin if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);} 358d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 359bfa639d9SMatthew G. Knepley PetscFunctionReturn(0); 360bfa639d9SMatthew G. Knepley } 361bfa639d9SMatthew G. Knepley 36289710940SMatthew G. Knepley /*@C 36389710940SMatthew G. Knepley PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 36489710940SMatthew G. Knepley 36589710940SMatthew G. Knepley Not collective 36689710940SMatthew G. Knepley 36789710940SMatthew G. Knepley Input Parameter: 36889710940SMatthew G. Knepley + q - The original PetscQuadrature 36989710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into 37089710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement 37189710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement 37289710940SMatthew G. Knepley 37389710940SMatthew G. Knepley Output Parameters: 37489710940SMatthew G. Knepley . dim - The dimension 37589710940SMatthew G. Knepley 37689710940SMatthew G. Knepley Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 37789710940SMatthew G. Knepley 378f5f57ec0SBarry Smith Not available from Fortran 379f5f57ec0SBarry Smith 38089710940SMatthew G. Knepley Level: intermediate 38189710940SMatthew G. Knepley 38289710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 38389710940SMatthew G. Knepley @*/ 38489710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 38589710940SMatthew G. Knepley { 38689710940SMatthew G. Knepley const PetscReal *points, *weights; 38789710940SMatthew G. Knepley PetscReal *pointsRef, *weightsRef; 388a6b92713SMatthew G. Knepley PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e; 38989710940SMatthew G. Knepley PetscErrorCode ierr; 39089710940SMatthew G. Knepley 39189710940SMatthew G. Knepley PetscFunctionBegin; 39289710940SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 39389710940SMatthew G. Knepley PetscValidPointer(v0, 3); 39489710940SMatthew G. Knepley PetscValidPointer(jac, 4); 39589710940SMatthew G. Knepley PetscValidPointer(qref, 5); 39689710940SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 39789710940SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 398a6b92713SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 39989710940SMatthew G. Knepley npointsRef = npoints*numSubelements; 40089710940SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 401a6b92713SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr); 40289710940SMatthew G. Knepley for (c = 0; c < numSubelements; ++c) { 40389710940SMatthew G. Knepley for (p = 0; p < npoints; ++p) { 40489710940SMatthew G. Knepley for (d = 0; d < dim; ++d) { 40589710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 40689710940SMatthew G. Knepley for (e = 0; e < dim; ++e) { 40789710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 40889710940SMatthew G. Knepley } 40989710940SMatthew G. Knepley } 41089710940SMatthew G. Knepley /* Could also use detJ here */ 411a6b92713SMatthew G. Knepley for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements; 41289710940SMatthew G. Knepley } 41389710940SMatthew G. Knepley } 41489710940SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 415a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 41689710940SMatthew G. Knepley PetscFunctionReturn(0); 41789710940SMatthew G. Knepley } 41889710940SMatthew G. Knepley 41937045ce4SJed Brown /*@ 42037045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 42137045ce4SJed Brown 42237045ce4SJed Brown Not Collective 42337045ce4SJed Brown 42437045ce4SJed Brown Input Arguments: 42537045ce4SJed Brown + npoints - number of spatial points to evaluate at 42637045ce4SJed Brown . points - array of locations to evaluate at 42737045ce4SJed Brown . ndegree - number of basis degrees to evaluate 42837045ce4SJed Brown - degrees - sorted array of degrees to evaluate 42937045ce4SJed Brown 43037045ce4SJed Brown Output Arguments: 4310298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 4320298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 4330298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 43437045ce4SJed Brown 43537045ce4SJed Brown Level: intermediate 43637045ce4SJed Brown 43737045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 43837045ce4SJed Brown @*/ 43937045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 44037045ce4SJed Brown { 44137045ce4SJed Brown PetscInt i,maxdegree; 44237045ce4SJed Brown 44337045ce4SJed Brown PetscFunctionBegin; 44437045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 44537045ce4SJed Brown maxdegree = degrees[ndegree-1]; 44637045ce4SJed Brown for (i=0; i<npoints; i++) { 44737045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 44837045ce4SJed Brown PetscInt j,k; 44937045ce4SJed Brown x = points[i]; 45037045ce4SJed Brown pm2 = 0; 45137045ce4SJed Brown pm1 = 1; 45237045ce4SJed Brown pd2 = 0; 45337045ce4SJed Brown pd1 = 0; 45437045ce4SJed Brown pdd2 = 0; 45537045ce4SJed Brown pdd1 = 0; 45637045ce4SJed Brown k = 0; 45737045ce4SJed Brown if (degrees[k] == 0) { 45837045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 45937045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 46037045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 46137045ce4SJed Brown k++; 46237045ce4SJed Brown } 46337045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 46437045ce4SJed Brown PetscReal p,d,dd; 46537045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 46637045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 46737045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 46837045ce4SJed Brown pm2 = pm1; 46937045ce4SJed Brown pm1 = p; 47037045ce4SJed Brown pd2 = pd1; 47137045ce4SJed Brown pd1 = d; 47237045ce4SJed Brown pdd2 = pdd1; 47337045ce4SJed Brown pdd1 = dd; 47437045ce4SJed Brown if (degrees[k] == j) { 47537045ce4SJed Brown if (B) B[i*ndegree+k] = p; 47637045ce4SJed Brown if (D) D[i*ndegree+k] = d; 47737045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 47837045ce4SJed Brown } 47937045ce4SJed Brown } 48037045ce4SJed Brown } 48137045ce4SJed Brown PetscFunctionReturn(0); 48237045ce4SJed Brown } 48337045ce4SJed Brown 48437045ce4SJed Brown /*@ 48537045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 48637045ce4SJed Brown 48737045ce4SJed Brown Not Collective 48837045ce4SJed Brown 48937045ce4SJed Brown Input Arguments: 49037045ce4SJed Brown + npoints - number of points 49137045ce4SJed Brown . a - left end of interval (often-1) 49237045ce4SJed Brown - b - right end of interval (often +1) 49337045ce4SJed Brown 49437045ce4SJed Brown Output Arguments: 49537045ce4SJed Brown + x - quadrature points 49637045ce4SJed Brown - w - quadrature weights 49737045ce4SJed Brown 49837045ce4SJed Brown Level: intermediate 49937045ce4SJed Brown 50037045ce4SJed Brown References: 50196a0c994SBarry Smith . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 50237045ce4SJed Brown 50337045ce4SJed Brown .seealso: PetscDTLegendreEval() 50437045ce4SJed Brown @*/ 50537045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 50637045ce4SJed Brown { 50737045ce4SJed Brown PetscErrorCode ierr; 50837045ce4SJed Brown PetscInt i; 50937045ce4SJed Brown PetscReal *work; 51037045ce4SJed Brown PetscScalar *Z; 51137045ce4SJed Brown PetscBLASInt N,LDZ,info; 51237045ce4SJed Brown 51337045ce4SJed Brown PetscFunctionBegin; 5140bfcf5a5SMatthew G. Knepley ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 51537045ce4SJed Brown /* Set up the Golub-Welsch system */ 51637045ce4SJed Brown for (i=0; i<npoints; i++) { 51737045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 51837045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 51937045ce4SJed Brown } 520dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 521c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 52237045ce4SJed Brown LDZ = N; 52337045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 5248b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 52537045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 5261c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 52737045ce4SJed Brown 52837045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 52937045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 53037045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 53119a57d60SBarry Smith if (x[i] == -0.0) x[i] = 0.0; 53237045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 5330d644c17SKarl Rupp 53488393a60SJed Brown w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 53537045ce4SJed Brown } 53637045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 53737045ce4SJed Brown PetscFunctionReturn(0); 53837045ce4SJed Brown } 539194825f6SJed Brown 5408272889dSSatish Balay static void qAndLEvaluation(PetscInt n, PetscReal x, PetscReal *q, PetscReal *qp, PetscReal *Ln) 5418272889dSSatish Balay /* 5428272889dSSatish Balay Compute the polynomial q(x) = L_{N+1}(x) - L_{n-1}(x) and its derivative in 5438272889dSSatish Balay addition to L_N(x) as these are needed for computing the GLL points via Newton's method. 5448272889dSSatish Balay Reference: "Implementing Spectral Methods for Partial Differential Equations: Algorithms 5458272889dSSatish Balay for Scientists and Engineers" by David A. Kopriva. 5468272889dSSatish Balay */ 5478272889dSSatish Balay { 5488272889dSSatish Balay PetscInt k; 5498272889dSSatish Balay 5508272889dSSatish Balay PetscReal Lnp; 5518272889dSSatish Balay PetscReal Lnp1, Lnp1p; 5528272889dSSatish Balay PetscReal Lnm1, Lnm1p; 5538272889dSSatish Balay PetscReal Lnm2, Lnm2p; 5548272889dSSatish Balay 5558272889dSSatish Balay Lnm1 = 1.0; 5568272889dSSatish Balay *Ln = x; 5578272889dSSatish Balay Lnm1p = 0.0; 5588272889dSSatish Balay Lnp = 1.0; 5598272889dSSatish Balay 5608272889dSSatish Balay for (k=2; k<=n; ++k) { 5618272889dSSatish Balay Lnm2 = Lnm1; 5628272889dSSatish Balay Lnm1 = *Ln; 5638272889dSSatish Balay Lnm2p = Lnm1p; 5648272889dSSatish Balay Lnm1p = Lnp; 5658272889dSSatish Balay *Ln = (2.*((PetscReal)k)-1.)/(1.0*((PetscReal)k))*x*Lnm1 - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm2; 5668272889dSSatish Balay Lnp = Lnm2p + (2.0*((PetscReal)k)-1.)*Lnm1; 5678272889dSSatish Balay } 5688272889dSSatish Balay k = n+1; 5698272889dSSatish Balay Lnp1 = (2.*((PetscReal)k)-1.)/(((PetscReal)k))*x*(*Ln) - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm1; 5708272889dSSatish Balay Lnp1p = Lnm1p + (2.0*((PetscReal)k)-1.)*(*Ln); 5718272889dSSatish Balay *q = Lnp1 - Lnm1; 5728272889dSSatish Balay *qp = Lnp1p - Lnm1p; 5738272889dSSatish Balay } 5748272889dSSatish Balay 5758272889dSSatish Balay /*@C 5768272889dSSatish Balay PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre 5778272889dSSatish Balay nodes of a given size on the domain [-1,1] 5788272889dSSatish Balay 5798272889dSSatish Balay Not Collective 5808272889dSSatish Balay 5818272889dSSatish Balay Input Parameter: 5828272889dSSatish Balay + n - number of grid nodes 583f2e8fe4dShannah_mairs - type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 5848272889dSSatish Balay 5858272889dSSatish Balay Output Arguments: 5868272889dSSatish Balay + x - quadrature points 5878272889dSSatish Balay - w - quadrature weights 5888272889dSSatish Balay 5898272889dSSatish Balay Notes: 5908272889dSSatish Balay For n > 30 the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not 5918272889dSSatish Balay close enough to the desired solution 5928272889dSSatish Balay 5938272889dSSatish Balay These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes 5948272889dSSatish Balay 595a8d69d7bSBarry 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 5968272889dSSatish Balay 5978272889dSSatish Balay Level: intermediate 5988272889dSSatish Balay 5998272889dSSatish Balay .seealso: PetscDTGaussQuadrature() 6008272889dSSatish Balay 6018272889dSSatish Balay @*/ 602916e780bShannah_mairs PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w) 6038272889dSSatish Balay { 6048272889dSSatish Balay PetscErrorCode ierr; 6058272889dSSatish Balay 6068272889dSSatish Balay PetscFunctionBegin; 6078272889dSSatish Balay if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element"); 6088272889dSSatish Balay 609f2e8fe4dShannah_mairs if (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA) { 6108272889dSSatish Balay PetscReal *M,si; 6118272889dSSatish Balay PetscBLASInt bn,lierr; 6128272889dSSatish Balay PetscReal x0,z0,z1,z2; 6138272889dSSatish Balay PetscInt i,p = npoints - 1,nn; 6148272889dSSatish Balay 6158272889dSSatish Balay x[0] =-1.0; 6168272889dSSatish Balay x[npoints-1] = 1.0; 6178272889dSSatish Balay if (npoints-2 > 0){ 6188272889dSSatish Balay ierr = PetscMalloc1(npoints-1,&M);CHKERRQ(ierr); 6198272889dSSatish Balay for (i=0; i<npoints-2; i++) { 6208272889dSSatish Balay si = ((PetscReal)i)+1.0; 6218272889dSSatish Balay M[i]=0.5*PetscSqrtReal(si*(si+2.0)/((si+0.5)*(si+1.5))); 6228272889dSSatish Balay } 6238272889dSSatish Balay ierr = PetscBLASIntCast(npoints-2,&bn);CHKERRQ(ierr); 6248272889dSSatish Balay ierr = PetscMemzero(&x[1],bn*sizeof(x[1]));CHKERRQ(ierr); 6258272889dSSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 6268272889dSSatish Balay x0=0; 6278272889dSSatish Balay PetscStackCallBLAS("LAPACKsteqr",LAPACKREALsteqr_("N",&bn,&x[1],M,&x0,&bn,M,&lierr)); 6288272889dSSatish Balay if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in STERF Lapack routine %d",(int)lierr); 6298272889dSSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 6308272889dSSatish Balay ierr = PetscFree(M);CHKERRQ(ierr); 6318272889dSSatish Balay } 6328272889dSSatish Balay if ((npoints-1)%2==0) { 6338272889dSSatish Balay x[(npoints-1)/2] = 0.0; /* hard wire to exactly 0.0 since linear algebra produces nonzero */ 6348272889dSSatish Balay } 6358272889dSSatish Balay 6368272889dSSatish Balay w[0] = w[p] = 2.0/(((PetscReal)(p))*(((PetscReal)p)+1.0)); 6378272889dSSatish Balay z2 = -1.; /* Dummy value to avoid -Wmaybe-initialized */ 6388272889dSSatish Balay for (i=1; i<p; i++) { 6398272889dSSatish Balay x0 = x[i]; 6408272889dSSatish Balay z0 = 1.0; 6418272889dSSatish Balay z1 = x0; 6428272889dSSatish Balay for (nn=1; nn<p; nn++) { 6438272889dSSatish Balay z2 = x0*z1*(2.0*((PetscReal)nn)+1.0)/(((PetscReal)nn)+1.0)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.0)); 6448272889dSSatish Balay z0 = z1; 6458272889dSSatish Balay z1 = z2; 6468272889dSSatish Balay } 6478272889dSSatish Balay w[i]=2.0/(((PetscReal)p)*(((PetscReal)p)+1.0)*z2*z2); 6488272889dSSatish Balay } 6498272889dSSatish Balay } else { 6508272889dSSatish Balay PetscInt j,m; 6518272889dSSatish Balay PetscReal z1,z,q,qp,Ln; 6528272889dSSatish Balay PetscReal *pt; 6538272889dSSatish Balay ierr = PetscMalloc1(npoints,&pt);CHKERRQ(ierr); 6548272889dSSatish Balay 655d410ae54Shannah_mairs if (npoints > 30) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON produces incorrect answers for n > 30"); 6568272889dSSatish Balay x[0] = -1.0; 6578272889dSSatish Balay x[npoints-1] = 1.0; 6588272889dSSatish Balay w[0] = w[npoints-1] = 2./(((PetscReal)npoints)*(((PetscReal)npoints)-1.0));; 6598272889dSSatish Balay m = (npoints-1)/2; /* The roots are symmetric, so we only find half of them. */ 6608272889dSSatish Balay for (j=1; j<=m; j++) { /* Loop over the desired roots. */ 6618272889dSSatish 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)); 6628272889dSSatish Balay /* Starting with the above approximation to the ith root, we enter */ 6638272889dSSatish Balay /* the main loop of refinement by Newton's method. */ 6648272889dSSatish Balay do { 6658272889dSSatish Balay qAndLEvaluation(npoints-1,z,&q,&qp,&Ln); 6668272889dSSatish Balay z1 = z; 6678272889dSSatish Balay z = z1-q/qp; /* Newton's method. */ 6688272889dSSatish Balay } while (PetscAbs(z-z1) > 10.*PETSC_MACHINE_EPSILON); 6698272889dSSatish Balay qAndLEvaluation(npoints-1,z,&q,&qp,&Ln); 6708272889dSSatish Balay 6718272889dSSatish Balay x[j] = z; 6728272889dSSatish Balay x[npoints-1-j] = -z; /* and put in its symmetric counterpart. */ 6738272889dSSatish Balay w[j] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln); /* Compute the weight */ 6748272889dSSatish Balay w[npoints-1-j] = w[j]; /* and its symmetric counterpart. */ 6758272889dSSatish Balay pt[j]=qp; 6768272889dSSatish Balay } 6778272889dSSatish Balay 6788272889dSSatish Balay if ((npoints-1)%2==0) { 6798272889dSSatish Balay qAndLEvaluation(npoints-1,0.0,&q,&qp,&Ln); 6808272889dSSatish Balay x[(npoints-1)/2] = 0.0; 6818272889dSSatish Balay w[(npoints-1)/2] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln); 6828272889dSSatish Balay } 6838272889dSSatish Balay ierr = PetscFree(pt);CHKERRQ(ierr); 6848272889dSSatish Balay } 6858272889dSSatish Balay PetscFunctionReturn(0); 6868272889dSSatish Balay } 6878272889dSSatish Balay 688744bafbcSMatthew G. Knepley /*@ 689744bafbcSMatthew G. Knepley PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 690744bafbcSMatthew G. Knepley 691744bafbcSMatthew G. Knepley Not Collective 692744bafbcSMatthew G. Knepley 693744bafbcSMatthew G. Knepley Input Arguments: 694744bafbcSMatthew G. Knepley + dim - The spatial dimension 695a6b92713SMatthew G. Knepley . Nc - The number of components 696744bafbcSMatthew G. Knepley . npoints - number of points in one dimension 697744bafbcSMatthew G. Knepley . a - left end of interval (often-1) 698744bafbcSMatthew G. Knepley - b - right end of interval (often +1) 699744bafbcSMatthew G. Knepley 700744bafbcSMatthew G. Knepley Output Argument: 701744bafbcSMatthew G. Knepley . q - A PetscQuadrature object 702744bafbcSMatthew G. Knepley 703744bafbcSMatthew G. Knepley Level: intermediate 704744bafbcSMatthew G. Knepley 705744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 706744bafbcSMatthew G. Knepley @*/ 707a6b92713SMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 708744bafbcSMatthew G. Knepley { 709a6b92713SMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c; 710744bafbcSMatthew G. Knepley PetscReal *x, *w, *xw, *ww; 711744bafbcSMatthew G. Knepley PetscErrorCode ierr; 712744bafbcSMatthew G. Knepley 713744bafbcSMatthew G. Knepley PetscFunctionBegin; 714744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 715a6b92713SMatthew G. Knepley ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr); 716744bafbcSMatthew G. Knepley /* Set up the Golub-Welsch system */ 717744bafbcSMatthew G. Knepley switch (dim) { 718744bafbcSMatthew G. Knepley case 0: 719744bafbcSMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 720744bafbcSMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 721744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 722a6b92713SMatthew G. Knepley ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 723744bafbcSMatthew G. Knepley x[0] = 0.0; 724a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[c] = 1.0; 725744bafbcSMatthew G. Knepley break; 726744bafbcSMatthew G. Knepley case 1: 727a6b92713SMatthew G. Knepley ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr); 728a6b92713SMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr); 729a6b92713SMatthew G. Knepley for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i]; 730a6b92713SMatthew G. Knepley ierr = PetscFree(ww);CHKERRQ(ierr); 731744bafbcSMatthew G. Knepley break; 732744bafbcSMatthew G. Knepley case 2: 733744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 734744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 735744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 736744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 737744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+0] = xw[i]; 738744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+1] = xw[j]; 739a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j]; 740744bafbcSMatthew G. Knepley } 741744bafbcSMatthew G. Knepley } 742744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 743744bafbcSMatthew G. Knepley break; 744744bafbcSMatthew G. Knepley case 3: 745744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 746744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 747744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 748744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 749744bafbcSMatthew G. Knepley for (k = 0; k < npoints; ++k) { 750744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 751744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 752744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 753a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k]; 754744bafbcSMatthew G. Knepley } 755744bafbcSMatthew G. Knepley } 756744bafbcSMatthew G. Knepley } 757744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 758744bafbcSMatthew G. Knepley break; 759744bafbcSMatthew G. Knepley default: 760744bafbcSMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 761744bafbcSMatthew G. Knepley } 762744bafbcSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 7632f5fb066SToby Isaac ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 764a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 765d9bac1caSLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr); 766744bafbcSMatthew G. Knepley PetscFunctionReturn(0); 767744bafbcSMatthew G. Knepley } 768744bafbcSMatthew G. Knepley 769494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 770494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 771494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 772494e7359SMatthew G. Knepley { 773494e7359SMatthew G. Knepley PetscReal f = 1.0; 774494e7359SMatthew G. Knepley PetscInt i; 775494e7359SMatthew G. Knepley 776494e7359SMatthew G. Knepley PetscFunctionBegin; 777494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 778494e7359SMatthew G. Knepley *factorial = f; 779494e7359SMatthew G. Knepley PetscFunctionReturn(0); 780494e7359SMatthew G. Knepley } 781494e7359SMatthew G. Knepley 782494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 783494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 784494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 785494e7359SMatthew G. Knepley { 786494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 787494e7359SMatthew G. Knepley PetscInt k; 788494e7359SMatthew G. Knepley 789494e7359SMatthew G. Knepley PetscFunctionBegin; 790494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 791494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 792494e7359SMatthew G. Knepley apb = a + b; 793494e7359SMatthew G. Knepley pn2 = 1.0; 794494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 795494e7359SMatthew G. Knepley *P = 0.0; 796494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 797494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 798494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 799494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 800494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 801494e7359SMatthew G. Knepley 802494e7359SMatthew G. Knepley a2 = a2 / a1; 803494e7359SMatthew G. Knepley a3 = a3 / a1; 804494e7359SMatthew G. Knepley a4 = a4 / a1; 805494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 806494e7359SMatthew G. Knepley pn2 = pn1; 807494e7359SMatthew G. Knepley pn1 = *P; 808494e7359SMatthew G. Knepley } 809494e7359SMatthew G. Knepley PetscFunctionReturn(0); 810494e7359SMatthew G. Knepley } 811494e7359SMatthew G. Knepley 812494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 813494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 814494e7359SMatthew G. Knepley { 815494e7359SMatthew G. Knepley PetscReal nP; 816494e7359SMatthew G. Knepley PetscErrorCode ierr; 817494e7359SMatthew G. Knepley 818494e7359SMatthew G. Knepley PetscFunctionBegin; 819494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 820494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 821494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 822494e7359SMatthew G. Knepley PetscFunctionReturn(0); 823494e7359SMatthew G. Knepley } 824494e7359SMatthew G. Knepley 825494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 826494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 827494e7359SMatthew G. Knepley { 828494e7359SMatthew G. Knepley PetscFunctionBegin; 829494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 830494e7359SMatthew G. Knepley *eta = y; 831494e7359SMatthew G. Knepley PetscFunctionReturn(0); 832494e7359SMatthew G. Knepley } 833494e7359SMatthew G. Knepley 834494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 835494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 836494e7359SMatthew G. Knepley { 837494e7359SMatthew G. Knepley PetscFunctionBegin; 838494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 839494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 840494e7359SMatthew G. Knepley *zeta = z; 841494e7359SMatthew G. Knepley PetscFunctionReturn(0); 842494e7359SMatthew G. Knepley } 843494e7359SMatthew G. Knepley 844494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 845494e7359SMatthew G. Knepley { 846494e7359SMatthew G. Knepley PetscInt maxIter = 100; 847494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 848a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 849494e7359SMatthew G. Knepley PetscInt k; 850494e7359SMatthew G. Knepley PetscErrorCode ierr; 851494e7359SMatthew G. Knepley 852494e7359SMatthew G. Knepley PetscFunctionBegin; 853a8291ba1SSatish Balay 8548b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 855a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 8560646a658SBarry Smith a2 = PetscTGamma(a + npoints + 1); 8570646a658SBarry Smith a3 = PetscTGamma(b + npoints + 1); 8580646a658SBarry Smith a4 = PetscTGamma(a + b + npoints + 1); 859a8291ba1SSatish Balay #else 86029bcbfd0SToby Isaac { 861d24bbb91SToby Isaac PetscInt ia, ib; 86229bcbfd0SToby Isaac 863d24bbb91SToby Isaac ia = (PetscInt) a; 864d24bbb91SToby Isaac ib = (PetscInt) b; 865d24bbb91SToby 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 */ 866d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr); 867d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr); 868d24bbb91SToby Isaac ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr); 86929bcbfd0SToby Isaac } else { 870a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 87129bcbfd0SToby Isaac } 87229bcbfd0SToby Isaac } 873a8291ba1SSatish Balay #endif 874a8291ba1SSatish Balay 875494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 876494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 877494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 878494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 879494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 8808b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 881494e7359SMatthew G. Knepley PetscInt j; 882494e7359SMatthew G. Knepley 883494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 884494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 885494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 886494e7359SMatthew G. Knepley PetscInt i; 887494e7359SMatthew G. Knepley 888494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 889494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 890494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 891494e7359SMatthew G. Knepley delta = f / (fp - f * s); 892494e7359SMatthew G. Knepley r = r - delta; 89377b4d14cSPeter Brune if (PetscAbsReal(delta) < eps) break; 894494e7359SMatthew G. Knepley } 895494e7359SMatthew G. Knepley x[k] = r; 896494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 897494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 898494e7359SMatthew G. Knepley } 899494e7359SMatthew G. Knepley PetscFunctionReturn(0); 900494e7359SMatthew G. Knepley } 901494e7359SMatthew G. Knepley 902f5f57ec0SBarry Smith /*@ 903494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 904494e7359SMatthew G. Knepley 905494e7359SMatthew G. Knepley Not Collective 906494e7359SMatthew G. Knepley 907494e7359SMatthew G. Knepley Input Arguments: 908494e7359SMatthew G. Knepley + dim - The simplex dimension 909a6b92713SMatthew G. Knepley . Nc - The number of components 910dcce0ee2SMatthew G. Knepley . npoints - The number of points in one dimension 911494e7359SMatthew G. Knepley . a - left end of interval (often-1) 912494e7359SMatthew G. Knepley - b - right end of interval (often +1) 913494e7359SMatthew G. Knepley 914744bafbcSMatthew G. Knepley Output Argument: 915552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 916494e7359SMatthew G. Knepley 917494e7359SMatthew G. Knepley Level: intermediate 918494e7359SMatthew G. Knepley 919494e7359SMatthew G. Knepley References: 92096a0c994SBarry Smith . 1. - Karniadakis and Sherwin. FIAT 921494e7359SMatthew G. Knepley 922744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 923494e7359SMatthew G. Knepley @*/ 924dcce0ee2SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 925494e7359SMatthew G. Knepley { 926dcce0ee2SMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints; 927494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 928a6b92713SMatthew G. Knepley PetscInt i, j, k, c; 929494e7359SMatthew G. Knepley PetscErrorCode ierr; 930494e7359SMatthew G. Knepley 931494e7359SMatthew G. Knepley PetscFunctionBegin; 932494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 933dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr); 934dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr); 935494e7359SMatthew G. Knepley switch (dim) { 936707aa5c5SMatthew G. Knepley case 0: 937707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 938707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 939785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 940a6b92713SMatthew G. Knepley ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 941707aa5c5SMatthew G. Knepley x[0] = 0.0; 942a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[c] = 1.0; 943707aa5c5SMatthew G. Knepley break; 944494e7359SMatthew G. Knepley case 1: 945dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr); 946dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr); 947dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i]; 948a6b92713SMatthew G. Knepley ierr = PetscFree(wx);CHKERRQ(ierr); 949494e7359SMatthew G. Knepley break; 950494e7359SMatthew G. Knepley case 2: 951dcce0ee2SMatthew G. Knepley ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr); 952dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 953dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 954dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 955dcce0ee2SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 956dcce0ee2SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 957dcce0ee2SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j]; 958494e7359SMatthew G. Knepley } 959494e7359SMatthew G. Knepley } 960494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 961494e7359SMatthew G. Knepley break; 962494e7359SMatthew G. Knepley case 3: 963dcce0ee2SMatthew G. Knepley ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr); 964dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 965dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 966dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 967dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 968dcce0ee2SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 969dcce0ee2SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 970dcce0ee2SMatthew 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); 971dcce0ee2SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k]; 972494e7359SMatthew G. Knepley } 973494e7359SMatthew G. Knepley } 974494e7359SMatthew G. Knepley } 975494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 976494e7359SMatthew G. Knepley break; 977494e7359SMatthew G. Knepley default: 978494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 979494e7359SMatthew G. Knepley } 98021454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 9812f5fb066SToby Isaac ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 982dcce0ee2SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 983d9bac1caSLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr); 984494e7359SMatthew G. Knepley PetscFunctionReturn(0); 985494e7359SMatthew G. Knepley } 986494e7359SMatthew G. Knepley 987f5f57ec0SBarry Smith /*@ 988b3c0f97bSTom Klotz PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 989b3c0f97bSTom Klotz 990b3c0f97bSTom Klotz Not Collective 991b3c0f97bSTom Klotz 992b3c0f97bSTom Klotz Input Arguments: 993b3c0f97bSTom Klotz + dim - The cell dimension 994b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l 995b3c0f97bSTom Klotz . a - left end of interval (often-1) 996b3c0f97bSTom Klotz - b - right end of interval (often +1) 997b3c0f97bSTom Klotz 998b3c0f97bSTom Klotz Output Argument: 999b3c0f97bSTom Klotz . q - A PetscQuadrature object 1000b3c0f97bSTom Klotz 1001b3c0f97bSTom Klotz Level: intermediate 1002b3c0f97bSTom Klotz 1003b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature() 1004b3c0f97bSTom Klotz @*/ 1005b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 1006b3c0f97bSTom Klotz { 1007b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 1008b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1009b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1010b3c0f97bSTom Klotz const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 1011d84b4d08SMatthew G. Knepley PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 1012b3c0f97bSTom Klotz PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 1013b3c0f97bSTom Klotz PetscReal *x, *w; 1014b3c0f97bSTom Klotz PetscInt K, k, npoints; 1015b3c0f97bSTom Klotz PetscErrorCode ierr; 1016b3c0f97bSTom Klotz 1017b3c0f97bSTom Klotz PetscFunctionBegin; 1018b3c0f97bSTom Klotz if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 1019b3c0f97bSTom Klotz if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 1020b3c0f97bSTom Klotz /* Find K such that the weights are < 32 digits of precision */ 1021b3c0f97bSTom Klotz for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 10229add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 1023b3c0f97bSTom Klotz } 1024b3c0f97bSTom Klotz ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1025b3c0f97bSTom Klotz ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 1026b3c0f97bSTom Klotz npoints = 2*K-1; 1027b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 1028b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 1029b3c0f97bSTom Klotz /* Center term */ 1030b3c0f97bSTom Klotz x[0] = beta; 1031b3c0f97bSTom Klotz w[0] = 0.5*alpha*PETSC_PI; 1032b3c0f97bSTom Klotz for (k = 1; k < K; ++k) { 10339add2064SThomas Klotz wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 10341118d4bcSLisandro Dalcin xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h)); 1035b3c0f97bSTom Klotz x[2*k-1] = -alpha*xk+beta; 1036b3c0f97bSTom Klotz w[2*k-1] = wk; 1037b3c0f97bSTom Klotz x[2*k+0] = alpha*xk+beta; 1038b3c0f97bSTom Klotz w[2*k+0] = wk; 1039b3c0f97bSTom Klotz } 1040a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr); 1041b3c0f97bSTom Klotz PetscFunctionReturn(0); 1042b3c0f97bSTom Klotz } 1043b3c0f97bSTom Klotz 1044b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1045b3c0f97bSTom Klotz { 1046b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 1047b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1048b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1049b3c0f97bSTom Klotz PetscReal h = 1.0; /* Step size, length between x_k */ 1050b3c0f97bSTom Klotz PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 1051b3c0f97bSTom Klotz PetscReal osum = 0.0; /* Integral on last level */ 1052b3c0f97bSTom Klotz PetscReal psum = 0.0; /* Integral on the level before the last level */ 1053b3c0f97bSTom Klotz PetscReal sum; /* Integral on current level */ 1054446c295cSMatthew G. Knepley PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 1055b3c0f97bSTom Klotz PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 1056b3c0f97bSTom Klotz PetscReal wk; /* Quadrature weight at x_k */ 1057b3c0f97bSTom Klotz PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 1058b3c0f97bSTom Klotz PetscInt d; /* Digits of precision in the integral */ 1059b3c0f97bSTom Klotz 1060b3c0f97bSTom Klotz PetscFunctionBegin; 1061b3c0f97bSTom Klotz if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 1062b3c0f97bSTom Klotz /* Center term */ 1063b3c0f97bSTom Klotz func(beta, &lval); 1064b3c0f97bSTom Klotz sum = 0.5*alpha*PETSC_PI*lval; 1065b3c0f97bSTom Klotz /* */ 1066b3c0f97bSTom Klotz do { 1067b3c0f97bSTom Klotz PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 1068b3c0f97bSTom Klotz PetscInt k = 1; 1069b3c0f97bSTom Klotz 1070b3c0f97bSTom Klotz ++l; 1071b3c0f97bSTom Klotz /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 1072b3c0f97bSTom Klotz /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 1073b3c0f97bSTom Klotz psum = osum; 1074b3c0f97bSTom Klotz osum = sum; 1075b3c0f97bSTom Klotz h *= 0.5; 1076b3c0f97bSTom Klotz sum *= 0.5; 1077b3c0f97bSTom Klotz do { 10789add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1079446c295cSMatthew G. Knepley yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1080446c295cSMatthew G. Knepley lx = -alpha*(1.0 - yk)+beta; 1081446c295cSMatthew G. Knepley rx = alpha*(1.0 - yk)+beta; 1082b3c0f97bSTom Klotz func(lx, &lval); 1083b3c0f97bSTom Klotz func(rx, &rval); 1084b3c0f97bSTom Klotz lterm = alpha*wk*lval; 1085b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 1086b3c0f97bSTom Klotz sum += lterm; 1087b3c0f97bSTom Klotz rterm = alpha*wk*rval; 1088b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 1089b3c0f97bSTom Klotz sum += rterm; 1090b3c0f97bSTom Klotz ++k; 1091b3c0f97bSTom Klotz /* Only need to evaluate every other point on refined levels */ 1092b3c0f97bSTom Klotz if (l != 1) ++k; 10939add2064SThomas Klotz } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 1094b3c0f97bSTom Klotz 1095b3c0f97bSTom Klotz d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 1096b3c0f97bSTom Klotz d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 1097b3c0f97bSTom Klotz d3 = PetscLog10Real(maxTerm) - p; 109809d48545SBarry Smith if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 109909d48545SBarry Smith else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 1100b3c0f97bSTom Klotz d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 11019add2064SThomas Klotz } while (d < digits && l < 12); 1102b3c0f97bSTom Klotz *sol = sum; 1103e510cb1fSThomas Klotz 1104b3c0f97bSTom Klotz PetscFunctionReturn(0); 1105b3c0f97bSTom Klotz } 1106b3c0f97bSTom Klotz 1107497880caSRichard Tran Mills #if defined(PETSC_HAVE_MPFR) 110829f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 110929f144ccSMatthew G. Knepley { 1110e510cb1fSThomas Klotz const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 111129f144ccSMatthew G. Knepley PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 111229f144ccSMatthew G. Knepley mpfr_t alpha; /* Half-width of the integration interval */ 111329f144ccSMatthew G. Knepley mpfr_t beta; /* Center of the integration interval */ 111429f144ccSMatthew G. Knepley mpfr_t h; /* Step size, length between x_k */ 111529f144ccSMatthew G. Knepley mpfr_t osum; /* Integral on last level */ 111629f144ccSMatthew G. Knepley mpfr_t psum; /* Integral on the level before the last level */ 111729f144ccSMatthew G. Knepley mpfr_t sum; /* Integral on current level */ 111829f144ccSMatthew G. Knepley mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 111929f144ccSMatthew G. Knepley mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 112029f144ccSMatthew G. Knepley mpfr_t wk; /* Quadrature weight at x_k */ 112129f144ccSMatthew G. Knepley PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 112229f144ccSMatthew G. Knepley PetscInt d; /* Digits of precision in the integral */ 112329f144ccSMatthew G. Knepley mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 112429f144ccSMatthew G. Knepley 112529f144ccSMatthew G. Knepley PetscFunctionBegin; 112629f144ccSMatthew G. Knepley if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 112729f144ccSMatthew G. Knepley /* Create high precision storage */ 1128c9f744b5SMatthew 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); 112929f144ccSMatthew G. Knepley /* Initialization */ 113029f144ccSMatthew G. Knepley mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 113129f144ccSMatthew G. Knepley mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 113229f144ccSMatthew G. Knepley mpfr_set_d(osum, 0.0, MPFR_RNDN); 113329f144ccSMatthew G. Knepley mpfr_set_d(psum, 0.0, MPFR_RNDN); 113429f144ccSMatthew G. Knepley mpfr_set_d(h, 1.0, MPFR_RNDN); 113529f144ccSMatthew G. Knepley mpfr_const_pi(pi2, MPFR_RNDN); 113629f144ccSMatthew G. Knepley mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 113729f144ccSMatthew G. Knepley /* Center term */ 113829f144ccSMatthew G. Knepley func(0.5*(b+a), &lval); 113929f144ccSMatthew G. Knepley mpfr_set(sum, pi2, MPFR_RNDN); 114029f144ccSMatthew G. Knepley mpfr_mul(sum, sum, alpha, MPFR_RNDN); 114129f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 114229f144ccSMatthew G. Knepley /* */ 114329f144ccSMatthew G. Knepley do { 114429f144ccSMatthew G. Knepley PetscReal d1, d2, d3, d4; 114529f144ccSMatthew G. Knepley PetscInt k = 1; 114629f144ccSMatthew G. Knepley 114729f144ccSMatthew G. Knepley ++l; 114829f144ccSMatthew G. Knepley mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 114929f144ccSMatthew G. Knepley /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 115029f144ccSMatthew G. Knepley /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 115129f144ccSMatthew G. Knepley mpfr_set(psum, osum, MPFR_RNDN); 115229f144ccSMatthew G. Knepley mpfr_set(osum, sum, MPFR_RNDN); 115329f144ccSMatthew G. Knepley mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 115429f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 115529f144ccSMatthew G. Knepley do { 115629f144ccSMatthew G. Knepley mpfr_set_si(kh, k, MPFR_RNDN); 115729f144ccSMatthew G. Knepley mpfr_mul(kh, kh, h, MPFR_RNDN); 115829f144ccSMatthew G. Knepley /* Weight */ 115929f144ccSMatthew G. Knepley mpfr_set(wk, h, MPFR_RNDN); 116029f144ccSMatthew G. Knepley mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 116129f144ccSMatthew G. Knepley mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 116229f144ccSMatthew G. Knepley mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 116329f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 116429f144ccSMatthew G. Knepley mpfr_sqr(tmp, tmp, MPFR_RNDN); 116529f144ccSMatthew G. Knepley mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 116629f144ccSMatthew G. Knepley mpfr_div(wk, wk, tmp, MPFR_RNDN); 116729f144ccSMatthew G. Knepley /* Abscissa */ 116829f144ccSMatthew G. Knepley mpfr_set_d(yk, 1.0, MPFR_RNDZ); 116929f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 117029f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 117129f144ccSMatthew G. Knepley mpfr_exp(tmp, msinh, MPFR_RNDN); 117229f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 117329f144ccSMatthew G. Knepley /* Quadrature points */ 117429f144ccSMatthew G. Knepley mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 117529f144ccSMatthew G. Knepley mpfr_mul(lx, lx, alpha, MPFR_RNDU); 117629f144ccSMatthew G. Knepley mpfr_add(lx, lx, beta, MPFR_RNDU); 117729f144ccSMatthew G. Knepley mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 117829f144ccSMatthew G. Knepley mpfr_mul(rx, rx, alpha, MPFR_RNDD); 117929f144ccSMatthew G. Knepley mpfr_add(rx, rx, beta, MPFR_RNDD); 118029f144ccSMatthew G. Knepley /* Evaluation */ 118129f144ccSMatthew G. Knepley func(mpfr_get_d(lx, MPFR_RNDU), &lval); 118229f144ccSMatthew G. Knepley func(mpfr_get_d(rx, MPFR_RNDD), &rval); 118329f144ccSMatthew G. Knepley /* Update */ 118429f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 118529f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 118629f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 118729f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 118829f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 118929f144ccSMatthew G. Knepley mpfr_set(curTerm, tmp, MPFR_RNDN); 119029f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 119129f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 119229f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 119329f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 119429f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 119529f144ccSMatthew G. Knepley mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 119629f144ccSMatthew G. Knepley ++k; 119729f144ccSMatthew G. Knepley /* Only need to evaluate every other point on refined levels */ 119829f144ccSMatthew G. Knepley if (l != 1) ++k; 119929f144ccSMatthew G. Knepley mpfr_log10(tmp, wk, MPFR_RNDN); 120029f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 1201c9f744b5SMatthew G. Knepley } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 120229f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, osum, MPFR_RNDN); 120329f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 120429f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 120529f144ccSMatthew G. Knepley d1 = mpfr_get_d(tmp, MPFR_RNDN); 120629f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, psum, MPFR_RNDN); 120729f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 120829f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 120929f144ccSMatthew G. Knepley d2 = mpfr_get_d(tmp, MPFR_RNDN); 121029f144ccSMatthew G. Knepley mpfr_log10(tmp, maxTerm, MPFR_RNDN); 1211c9f744b5SMatthew G. Knepley d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 121229f144ccSMatthew G. Knepley mpfr_log10(tmp, curTerm, MPFR_RNDN); 121329f144ccSMatthew G. Knepley d4 = mpfr_get_d(tmp, MPFR_RNDN); 121429f144ccSMatthew G. Knepley d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1215b0649871SThomas Klotz } while (d < digits && l < 8); 121629f144ccSMatthew G. Knepley *sol = mpfr_get_d(sum, MPFR_RNDN); 121729f144ccSMatthew G. Knepley /* Cleanup */ 121829f144ccSMatthew G. Knepley mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 121929f144ccSMatthew G. Knepley PetscFunctionReturn(0); 122029f144ccSMatthew G. Knepley } 1221d525116cSMatthew G. Knepley #else 1222fbfcfee5SBarry Smith 1223d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1224d525116cSMatthew G. Knepley { 1225d525116cSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 1226d525116cSMatthew G. Knepley } 122729f144ccSMatthew G. Knepley #endif 122829f144ccSMatthew G. Knepley 1229194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 1230194825f6SJed Brown * A in column-major format 1231194825f6SJed Brown * Ainv in row-major format 1232194825f6SJed Brown * tau has length m 1233194825f6SJed Brown * worksize must be >= max(1,n) 1234194825f6SJed Brown */ 1235194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1236194825f6SJed Brown { 1237194825f6SJed Brown PetscErrorCode ierr; 1238194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1239194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 1240194825f6SJed Brown 1241194825f6SJed Brown PetscFunctionBegin; 1242194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1243194825f6SJed Brown { 1244194825f6SJed Brown PetscInt i,j; 1245dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 1246194825f6SJed Brown for (j=0; j<n; j++) { 1247194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 1248194825f6SJed Brown } 1249194825f6SJed Brown mstride = m; 1250194825f6SJed Brown } 1251194825f6SJed Brown #else 1252194825f6SJed Brown A = A_in; 1253194825f6SJed Brown Ainv = Ainv_out; 1254194825f6SJed Brown #endif 1255194825f6SJed Brown 1256194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 1257194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 1258194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 1259194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 1260194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1261001a771dSBarry Smith PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1262194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 1263194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1264194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1265194825f6SJed Brown 1266194825f6SJed Brown /* Extract an explicit representation of Q */ 1267194825f6SJed Brown Q = Ainv; 1268194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 1269194825f6SJed Brown K = N; /* full rank */ 1270c964aadfSJose E. Roman PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1271194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1272194825f6SJed Brown 1273194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1274194825f6SJed Brown Alpha = 1.0; 1275194825f6SJed Brown ldb = lda; 1276001a771dSBarry Smith PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 1277194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 1278194825f6SJed Brown 1279194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1280194825f6SJed Brown { 1281194825f6SJed Brown PetscInt i; 1282194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 1283194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 1284194825f6SJed Brown } 1285194825f6SJed Brown #endif 1286194825f6SJed Brown PetscFunctionReturn(0); 1287194825f6SJed Brown } 1288194825f6SJed Brown 1289194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 1290194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 1291194825f6SJed Brown { 1292194825f6SJed Brown PetscErrorCode ierr; 1293194825f6SJed Brown PetscReal *Bv; 1294194825f6SJed Brown PetscInt i,j; 1295194825f6SJed Brown 1296194825f6SJed Brown PetscFunctionBegin; 1297785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 1298194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 1299194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 1300194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 1301194825f6SJed Brown for (i=0; i<ninterval; i++) { 1302194825f6SJed Brown for (j=0; j<ndegree; j++) { 1303194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1304194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1305194825f6SJed Brown } 1306194825f6SJed Brown } 1307194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 1308194825f6SJed Brown PetscFunctionReturn(0); 1309194825f6SJed Brown } 1310194825f6SJed Brown 1311194825f6SJed Brown /*@ 1312194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 1313194825f6SJed Brown 1314194825f6SJed Brown Not Collective 1315194825f6SJed Brown 1316194825f6SJed Brown Input Arguments: 1317194825f6SJed Brown + degree - degree of reconstruction polynomial 1318194825f6SJed Brown . nsource - number of source intervals 1319194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 1320194825f6SJed Brown . ntarget - number of target intervals 1321194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 1322194825f6SJed Brown 1323194825f6SJed Brown Output Arguments: 1324194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 1325194825f6SJed Brown 1326194825f6SJed Brown Level: advanced 1327194825f6SJed Brown 1328194825f6SJed Brown .seealso: PetscDTLegendreEval() 1329194825f6SJed Brown @*/ 1330194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 1331194825f6SJed Brown { 1332194825f6SJed Brown PetscErrorCode ierr; 1333194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 1334194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 1335194825f6SJed Brown PetscScalar *tau,*work; 1336194825f6SJed Brown 1337194825f6SJed Brown PetscFunctionBegin; 1338194825f6SJed Brown PetscValidRealPointer(sourcex,3); 1339194825f6SJed Brown PetscValidRealPointer(targetx,5); 1340194825f6SJed Brown PetscValidRealPointer(R,6); 1341194825f6SJed 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); 1342194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 1343194825f6SJed Brown for (i=0; i<nsource; i++) { 134457622a8eSBarry 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]); 1345194825f6SJed Brown } 1346194825f6SJed Brown for (i=0; i<ntarget; i++) { 134757622a8eSBarry 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]); 1348194825f6SJed Brown } 1349194825f6SJed Brown #endif 1350194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 1351194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 1352194825f6SJed Brown center = (xmin + xmax)/2; 1353194825f6SJed Brown hscale = (xmax - xmin)/2; 1354194825f6SJed Brown worksize = nsource; 1355dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 1356dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 1357194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 1358194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 1359194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 1360194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 1361194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 1362194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 1363194825f6SJed Brown for (i=0; i<ntarget; i++) { 1364194825f6SJed Brown PetscReal rowsum = 0; 1365194825f6SJed Brown for (j=0; j<nsource; j++) { 1366194825f6SJed Brown PetscReal sum = 0; 1367194825f6SJed Brown for (k=0; k<degree+1; k++) { 1368194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 1369194825f6SJed Brown } 1370194825f6SJed Brown R[i*nsource+j] = sum; 1371194825f6SJed Brown rowsum += sum; 1372194825f6SJed Brown } 1373194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 1374194825f6SJed Brown } 1375194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 1376194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 1377194825f6SJed Brown PetscFunctionReturn(0); 1378194825f6SJed Brown } 1379916e780bShannah_mairs 1380916e780bShannah_mairs /*@C 1381916e780bShannah_mairs PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points 1382916e780bShannah_mairs 1383916e780bShannah_mairs Not Collective 1384916e780bShannah_mairs 1385916e780bShannah_mairs Input Parameter: 1386916e780bShannah_mairs + n - the number of GLL nodes 1387916e780bShannah_mairs . nodes - the GLL nodes 1388916e780bShannah_mairs . weights - the GLL weights 1389916e780bShannah_mairs . f - the function values at the nodes 1390916e780bShannah_mairs 1391916e780bShannah_mairs Output Parameter: 1392916e780bShannah_mairs . in - the value of the integral 1393916e780bShannah_mairs 1394916e780bShannah_mairs Level: beginner 1395916e780bShannah_mairs 1396916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature() 1397916e780bShannah_mairs 1398916e780bShannah_mairs @*/ 1399916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in) 1400916e780bShannah_mairs { 1401916e780bShannah_mairs PetscInt i; 1402916e780bShannah_mairs 1403916e780bShannah_mairs PetscFunctionBegin; 1404916e780bShannah_mairs *in = 0.; 1405916e780bShannah_mairs for (i=0; i<n; i++) { 1406916e780bShannah_mairs *in += f[i]*f[i]*weights[i]; 1407916e780bShannah_mairs } 1408916e780bShannah_mairs PetscFunctionReturn(0); 1409916e780bShannah_mairs } 1410916e780bShannah_mairs 1411916e780bShannah_mairs /*@C 1412916e780bShannah_mairs PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element 1413916e780bShannah_mairs 1414916e780bShannah_mairs Not Collective 1415916e780bShannah_mairs 1416916e780bShannah_mairs Input Parameter: 1417916e780bShannah_mairs + n - the number of GLL nodes 1418916e780bShannah_mairs . nodes - the GLL nodes 1419916e780bShannah_mairs . weights - the GLL weights 1420916e780bShannah_mairs 1421916e780bShannah_mairs Output Parameter: 1422916e780bShannah_mairs . A - the stiffness element 1423916e780bShannah_mairs 1424916e780bShannah_mairs Level: beginner 1425916e780bShannah_mairs 1426916e780bShannah_mairs Notes: 1427916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy() 1428916e780bShannah_mairs 1429916e780bShannah_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) 1430916e780bShannah_mairs 1431916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1432916e780bShannah_mairs 1433916e780bShannah_mairs @*/ 1434916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1435916e780bShannah_mairs { 1436916e780bShannah_mairs PetscReal **A; 1437916e780bShannah_mairs PetscErrorCode ierr; 1438916e780bShannah_mairs const PetscReal *gllnodes = nodes; 1439916e780bShannah_mairs const PetscInt p = n-1; 1440916e780bShannah_mairs PetscReal z0,z1,z2 = -1,x,Lpj,Lpr; 1441916e780bShannah_mairs PetscInt i,j,nn,r; 1442916e780bShannah_mairs 1443916e780bShannah_mairs PetscFunctionBegin; 1444916e780bShannah_mairs ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1445916e780bShannah_mairs ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1446916e780bShannah_mairs for (i=1; i<n; i++) A[i] = A[i-1]+n; 1447916e780bShannah_mairs 1448916e780bShannah_mairs for (j=1; j<p; j++) { 1449916e780bShannah_mairs x = gllnodes[j]; 1450916e780bShannah_mairs z0 = 1.; 1451916e780bShannah_mairs z1 = x; 1452916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1453916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1454916e780bShannah_mairs z0 = z1; 1455916e780bShannah_mairs z1 = z2; 1456916e780bShannah_mairs } 1457916e780bShannah_mairs Lpj=z2; 1458916e780bShannah_mairs for (r=1; r<p; r++) { 1459916e780bShannah_mairs if (r == j) { 1460916e780bShannah_mairs A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj); 1461916e780bShannah_mairs } else { 1462916e780bShannah_mairs x = gllnodes[r]; 1463916e780bShannah_mairs z0 = 1.; 1464916e780bShannah_mairs z1 = x; 1465916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1466916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1467916e780bShannah_mairs z0 = z1; 1468916e780bShannah_mairs z1 = z2; 1469916e780bShannah_mairs } 1470916e780bShannah_mairs Lpr = z2; 1471916e780bShannah_mairs A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r])); 1472916e780bShannah_mairs } 1473916e780bShannah_mairs } 1474916e780bShannah_mairs } 1475916e780bShannah_mairs for (j=1; j<p+1; j++) { 1476916e780bShannah_mairs x = gllnodes[j]; 1477916e780bShannah_mairs z0 = 1.; 1478916e780bShannah_mairs z1 = x; 1479916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1480916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1481916e780bShannah_mairs z0 = z1; 1482916e780bShannah_mairs z1 = z2; 1483916e780bShannah_mairs } 1484916e780bShannah_mairs Lpj = z2; 1485916e780bShannah_mairs A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j])); 1486916e780bShannah_mairs A[0][j] = A[j][0]; 1487916e780bShannah_mairs } 1488916e780bShannah_mairs for (j=0; j<p; j++) { 1489916e780bShannah_mairs x = gllnodes[j]; 1490916e780bShannah_mairs z0 = 1.; 1491916e780bShannah_mairs z1 = x; 1492916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1493916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1494916e780bShannah_mairs z0 = z1; 1495916e780bShannah_mairs z1 = z2; 1496916e780bShannah_mairs } 1497916e780bShannah_mairs Lpj=z2; 1498916e780bShannah_mairs 1499916e780bShannah_mairs A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j])); 1500916e780bShannah_mairs A[j][p] = A[p][j]; 1501916e780bShannah_mairs } 1502916e780bShannah_mairs A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.; 1503916e780bShannah_mairs A[p][p]=A[0][0]; 1504916e780bShannah_mairs *AA = A; 1505916e780bShannah_mairs PetscFunctionReturn(0); 1506916e780bShannah_mairs } 1507916e780bShannah_mairs 1508916e780bShannah_mairs /*@C 1509916e780bShannah_mairs PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element 1510916e780bShannah_mairs 1511916e780bShannah_mairs Not Collective 1512916e780bShannah_mairs 1513916e780bShannah_mairs Input Parameter: 1514916e780bShannah_mairs + n - the number of GLL nodes 1515916e780bShannah_mairs . nodes - the GLL nodes 1516916e780bShannah_mairs . weights - the GLL weightss 1517916e780bShannah_mairs - A - the stiffness element 1518916e780bShannah_mairs 1519916e780bShannah_mairs Level: beginner 1520916e780bShannah_mairs 1521916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate() 1522916e780bShannah_mairs 1523916e780bShannah_mairs @*/ 1524916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1525916e780bShannah_mairs { 1526916e780bShannah_mairs PetscErrorCode ierr; 1527916e780bShannah_mairs 1528916e780bShannah_mairs PetscFunctionBegin; 1529916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1530916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1531916e780bShannah_mairs *AA = NULL; 1532916e780bShannah_mairs PetscFunctionReturn(0); 1533916e780bShannah_mairs } 1534916e780bShannah_mairs 1535916e780bShannah_mairs /*@C 1536916e780bShannah_mairs PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element 1537916e780bShannah_mairs 1538916e780bShannah_mairs Not Collective 1539916e780bShannah_mairs 1540916e780bShannah_mairs Input Parameter: 1541916e780bShannah_mairs + n - the number of GLL nodes 1542916e780bShannah_mairs . nodes - the GLL nodes 1543916e780bShannah_mairs . weights - the GLL weights 1544916e780bShannah_mairs 1545916e780bShannah_mairs Output Parameter: 1546916e780bShannah_mairs . AA - the stiffness element 1547916e780bShannah_mairs - AAT - the transpose of AA (pass in NULL if you do not need this array) 1548916e780bShannah_mairs 1549916e780bShannah_mairs Level: beginner 1550916e780bShannah_mairs 1551916e780bShannah_mairs Notes: 1552916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementGradientDestroy() 1553916e780bShannah_mairs 1554916e780bShannah_mairs You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1555916e780bShannah_mairs 1556916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1557916e780bShannah_mairs 1558916e780bShannah_mairs @*/ 1559916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1560916e780bShannah_mairs { 1561916e780bShannah_mairs PetscReal **A, **AT = NULL; 1562916e780bShannah_mairs PetscErrorCode ierr; 1563916e780bShannah_mairs const PetscReal *gllnodes = nodes; 1564916e780bShannah_mairs const PetscInt p = n-1; 1565916e780bShannah_mairs PetscReal q,qp,Li, Lj,d0; 1566916e780bShannah_mairs PetscInt i,j; 1567916e780bShannah_mairs 1568916e780bShannah_mairs PetscFunctionBegin; 1569916e780bShannah_mairs ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1570916e780bShannah_mairs ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1571916e780bShannah_mairs for (i=1; i<n; i++) A[i] = A[i-1]+n; 1572916e780bShannah_mairs 1573916e780bShannah_mairs if (AAT) { 1574916e780bShannah_mairs ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr); 1575916e780bShannah_mairs ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr); 1576916e780bShannah_mairs for (i=1; i<n; i++) AT[i] = AT[i-1]+n; 1577916e780bShannah_mairs } 1578916e780bShannah_mairs 1579916e780bShannah_mairs if (n==1) {A[0][0] = 0.;} 1580916e780bShannah_mairs d0 = (PetscReal)p*((PetscReal)p+1.)/4.; 1581916e780bShannah_mairs for (i=0; i<n; i++) { 1582916e780bShannah_mairs for (j=0; j<n; j++) { 1583916e780bShannah_mairs A[i][j] = 0.; 1584fdd31e58Shannah_mairs qAndLEvaluation(p,gllnodes[i],&q,&qp,&Li); 1585fdd31e58Shannah_mairs qAndLEvaluation(p,gllnodes[j],&q,&qp,&Lj); 1586916e780bShannah_mairs if (i!=j) A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j])); 1587916e780bShannah_mairs if ((j==i) && (i==0)) A[i][j] = -d0; 1588916e780bShannah_mairs if (j==i && i==p) A[i][j] = d0; 1589916e780bShannah_mairs if (AT) AT[j][i] = A[i][j]; 1590916e780bShannah_mairs } 1591916e780bShannah_mairs } 1592916e780bShannah_mairs if (AAT) *AAT = AT; 1593916e780bShannah_mairs *AA = A; 1594916e780bShannah_mairs PetscFunctionReturn(0); 1595916e780bShannah_mairs } 1596916e780bShannah_mairs 1597916e780bShannah_mairs /*@C 1598916e780bShannah_mairs PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate() 1599916e780bShannah_mairs 1600916e780bShannah_mairs Not Collective 1601916e780bShannah_mairs 1602916e780bShannah_mairs Input Parameter: 1603916e780bShannah_mairs + n - the number of GLL nodes 1604916e780bShannah_mairs . nodes - the GLL nodes 1605916e780bShannah_mairs . weights - the GLL weights 1606916e780bShannah_mairs . AA - the stiffness element 1607916e780bShannah_mairs - AAT - the transpose of the element 1608916e780bShannah_mairs 1609916e780bShannah_mairs Level: beginner 1610916e780bShannah_mairs 1611916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate() 1612916e780bShannah_mairs 1613916e780bShannah_mairs @*/ 1614916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1615916e780bShannah_mairs { 1616916e780bShannah_mairs PetscErrorCode ierr; 1617916e780bShannah_mairs 1618916e780bShannah_mairs PetscFunctionBegin; 1619916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1620916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1621916e780bShannah_mairs *AA = NULL; 1622916e780bShannah_mairs if (*AAT) { 1623916e780bShannah_mairs ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr); 1624916e780bShannah_mairs ierr = PetscFree(*AAT);CHKERRQ(ierr); 1625916e780bShannah_mairs *AAT = NULL; 1626916e780bShannah_mairs } 1627916e780bShannah_mairs PetscFunctionReturn(0); 1628916e780bShannah_mairs } 1629916e780bShannah_mairs 1630916e780bShannah_mairs /*@C 1631916e780bShannah_mairs PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element 1632916e780bShannah_mairs 1633916e780bShannah_mairs Not Collective 1634916e780bShannah_mairs 1635916e780bShannah_mairs Input Parameter: 1636916e780bShannah_mairs + n - the number of GLL nodes 1637916e780bShannah_mairs . nodes - the GLL nodes 1638916e780bShannah_mairs . weights - the GLL weightss 1639916e780bShannah_mairs 1640916e780bShannah_mairs Output Parameter: 1641916e780bShannah_mairs . AA - the stiffness element 1642916e780bShannah_mairs 1643916e780bShannah_mairs Level: beginner 1644916e780bShannah_mairs 1645916e780bShannah_mairs Notes: 1646916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy() 1647916e780bShannah_mairs 1648916e780bShannah_mairs This is the same as the Gradient operator multiplied by the diagonal mass matrix 1649916e780bShannah_mairs 1650916e780bShannah_mairs You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1651916e780bShannah_mairs 1652916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy() 1653916e780bShannah_mairs 1654916e780bShannah_mairs @*/ 1655916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1656916e780bShannah_mairs { 1657916e780bShannah_mairs PetscReal **D; 1658916e780bShannah_mairs PetscErrorCode ierr; 1659916e780bShannah_mairs const PetscReal *gllweights = weights; 1660916e780bShannah_mairs const PetscInt glln = n; 1661916e780bShannah_mairs PetscInt i,j; 1662916e780bShannah_mairs 1663916e780bShannah_mairs PetscFunctionBegin; 1664916e780bShannah_mairs ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr); 1665916e780bShannah_mairs for (i=0; i<glln; i++){ 1666916e780bShannah_mairs for (j=0; j<glln; j++) { 1667916e780bShannah_mairs D[i][j] = gllweights[i]*D[i][j]; 1668916e780bShannah_mairs } 1669916e780bShannah_mairs } 1670916e780bShannah_mairs *AA = D; 1671916e780bShannah_mairs PetscFunctionReturn(0); 1672916e780bShannah_mairs } 1673916e780bShannah_mairs 1674916e780bShannah_mairs /*@C 1675916e780bShannah_mairs PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element 1676916e780bShannah_mairs 1677916e780bShannah_mairs Not Collective 1678916e780bShannah_mairs 1679916e780bShannah_mairs Input Parameter: 1680916e780bShannah_mairs + n - the number of GLL nodes 1681916e780bShannah_mairs . nodes - the GLL nodes 1682916e780bShannah_mairs . weights - the GLL weights 1683916e780bShannah_mairs - A - advection 1684916e780bShannah_mairs 1685916e780bShannah_mairs Level: beginner 1686916e780bShannah_mairs 1687916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate() 1688916e780bShannah_mairs 1689916e780bShannah_mairs @*/ 1690916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1691916e780bShannah_mairs { 1692916e780bShannah_mairs PetscErrorCode ierr; 1693916e780bShannah_mairs 1694916e780bShannah_mairs PetscFunctionBegin; 1695916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1696916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1697916e780bShannah_mairs *AA = NULL; 1698916e780bShannah_mairs PetscFunctionReturn(0); 1699916e780bShannah_mairs } 1700916e780bShannah_mairs 1701916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1702916e780bShannah_mairs { 1703916e780bShannah_mairs PetscReal **A; 1704916e780bShannah_mairs PetscErrorCode ierr; 1705916e780bShannah_mairs const PetscReal *gllweights = weights; 1706916e780bShannah_mairs const PetscInt glln = n; 1707916e780bShannah_mairs PetscInt i,j; 1708916e780bShannah_mairs 1709916e780bShannah_mairs PetscFunctionBegin; 1710916e780bShannah_mairs ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr); 1711916e780bShannah_mairs ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr); 1712916e780bShannah_mairs for (i=1; i<glln; i++) A[i] = A[i-1]+glln; 1713916e780bShannah_mairs if (glln==1) {A[0][0] = 0.;} 1714916e780bShannah_mairs for (i=0; i<glln; i++) { 1715916e780bShannah_mairs for (j=0; j<glln; j++) { 1716916e780bShannah_mairs A[i][j] = 0.; 1717916e780bShannah_mairs if (j==i) A[i][j] = gllweights[i]; 1718916e780bShannah_mairs } 1719916e780bShannah_mairs } 1720916e780bShannah_mairs *AA = A; 1721916e780bShannah_mairs PetscFunctionReturn(0); 1722916e780bShannah_mairs } 1723916e780bShannah_mairs 1724916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1725916e780bShannah_mairs { 1726916e780bShannah_mairs PetscErrorCode ierr; 1727916e780bShannah_mairs 1728916e780bShannah_mairs PetscFunctionBegin; 1729916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1730916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1731916e780bShannah_mairs *AA = NULL; 1732916e780bShannah_mairs PetscFunctionReturn(0); 1733916e780bShannah_mairs } 1734916e780bShannah_mairs 1735