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 28d083f849SBarry 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 60d083f849SBarry 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); 87580bdb30SBarry Smith ierr = PetscArraycpy(p, points, Nq*dim);CHKERRQ(ierr); 88580bdb30SBarry Smith ierr = PetscArraycpy(w, weights, Nc * Nq);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 96d083f849SBarry 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 265*907761f8SToby Isaac static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[]) 266*907761f8SToby Isaac { 267*907761f8SToby Isaac PetscScalar *Js, *Jinvs; 268*907761f8SToby Isaac PetscInt i, j, k; 269*907761f8SToby Isaac PetscBLASInt bm, bn, info; 270*907761f8SToby Isaac PetscErrorCode ierr; 271*907761f8SToby Isaac 272*907761f8SToby Isaac PetscFunctionBegin; 273*907761f8SToby Isaac ierr = PetscBLASIntCast(m, &bm);CHKERRQ(ierr); 274*907761f8SToby Isaac ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr); 275*907761f8SToby Isaac #if defined(PETSC_USE_COMPLEX) 276*907761f8SToby Isaac ierr = PetscMalloc2(m*n, &Js, m*n, &Jinvs);CHKERRQ(ierr); 277*907761f8SToby Isaac for (i = 0; i < m*n; j++) Js[i] = J[i]; 278*907761f8SToby Isaac #else 279*907761f8SToby Isaac Js = (PetscReal *) J; 280*907761f8SToby Isaac Jinvs = Jinv; 281*907761f8SToby Isaac #endif 282*907761f8SToby Isaac if (m == n) { 283*907761f8SToby Isaac PetscBLASInt *pivots; 284*907761f8SToby Isaac PetscScalar *W; 285*907761f8SToby Isaac 286*907761f8SToby Isaac ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr); 287*907761f8SToby Isaac 288*907761f8SToby Isaac ierr = PetscArraycpy(Jinvs, Js, m * m);CHKERRQ(ierr); 289*907761f8SToby Isaac PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info)); 290*907761f8SToby Isaac if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 291*907761f8SToby Isaac PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info)); 292*907761f8SToby Isaac if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 293*907761f8SToby Isaac ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 294*907761f8SToby Isaac } else if (m < n) { 295*907761f8SToby Isaac PetscScalar *JJT; 296*907761f8SToby Isaac PetscBLASInt *pivots; 297*907761f8SToby Isaac PetscScalar *W; 298*907761f8SToby Isaac 299*907761f8SToby Isaac ierr = PetscMalloc1(m*m, &JJT);CHKERRQ(ierr); 300*907761f8SToby Isaac ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr); 301*907761f8SToby Isaac for (i = 0; i < m; i++) { 302*907761f8SToby Isaac for (j = 0; j < m; j++) { 303*907761f8SToby Isaac PetscScalar val = 0.; 304*907761f8SToby Isaac 305*907761f8SToby Isaac for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k]; 306*907761f8SToby Isaac JJT[i * m + j] = val; 307*907761f8SToby Isaac } 308*907761f8SToby Isaac } 309*907761f8SToby Isaac 310*907761f8SToby Isaac PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info)); 311*907761f8SToby Isaac if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 312*907761f8SToby Isaac PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info)); 313*907761f8SToby Isaac if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 314*907761f8SToby Isaac for (i = 0; i < n; i++) { 315*907761f8SToby Isaac for (j = 0; j < m; j++) { 316*907761f8SToby Isaac PetscScalar val = 0.; 317*907761f8SToby Isaac 318*907761f8SToby Isaac for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j]; 319*907761f8SToby Isaac Jinvs[i * m + j] = val; 320*907761f8SToby Isaac } 321*907761f8SToby Isaac } 322*907761f8SToby Isaac ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 323*907761f8SToby Isaac ierr = PetscFree(JJT);CHKERRQ(ierr); 324*907761f8SToby Isaac } else { 325*907761f8SToby Isaac PetscScalar *JTJ; 326*907761f8SToby Isaac PetscBLASInt *pivots; 327*907761f8SToby Isaac PetscScalar *W; 328*907761f8SToby Isaac 329*907761f8SToby Isaac ierr = PetscMalloc1(n*n, &JTJ);CHKERRQ(ierr); 330*907761f8SToby Isaac ierr = PetscMalloc2(n, &pivots, n, &W);CHKERRQ(ierr); 331*907761f8SToby Isaac for (i = 0; i < n; i++) { 332*907761f8SToby Isaac for (j = 0; j < n; j++) { 333*907761f8SToby Isaac PetscScalar val = 0.; 334*907761f8SToby Isaac 335*907761f8SToby Isaac for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j]; 336*907761f8SToby Isaac JTJ[i * n + j] = val; 337*907761f8SToby Isaac } 338*907761f8SToby Isaac } 339*907761f8SToby Isaac 340*907761f8SToby Isaac PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bm, pivots, &info)); 341*907761f8SToby Isaac if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 342*907761f8SToby Isaac PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info)); 343*907761f8SToby Isaac if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 344*907761f8SToby Isaac for (i = 0; i < n; i++) { 345*907761f8SToby Isaac for (j = 0; j < m; j++) { 346*907761f8SToby Isaac PetscScalar val = 0.; 347*907761f8SToby Isaac 348*907761f8SToby Isaac for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k]; 349*907761f8SToby Isaac Jinvs[i * m + j] = val; 350*907761f8SToby Isaac } 351*907761f8SToby Isaac } 352*907761f8SToby Isaac ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 353*907761f8SToby Isaac ierr = PetscFree(JTJ);CHKERRQ(ierr); 354*907761f8SToby Isaac } 355*907761f8SToby Isaac #if defined(PETSC_USE_COMPLEX) 356*907761f8SToby Isaac for (i = 0; i < m*n; j++) Jinv[i] = PetscRealPart(Jinvs[i]); 357*907761f8SToby Isaac ierr = PetscFree2(Js, Jinvs);CHKERRQ(ierr); 358*907761f8SToby Isaac #endif 359*907761f8SToby Isaac PetscFunctionReturn(0); 360*907761f8SToby Isaac } 361*907761f8SToby Isaac 362*907761f8SToby Isaac /*@ 363*907761f8SToby Isaac PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation. 364*907761f8SToby Isaac 365*907761f8SToby Isaac Collecive on PetscQuadrature 366*907761f8SToby Isaac 367*907761f8SToby Isaac Input Arguments: 368*907761f8SToby Isaac + q - the quadrature functional 369*907761f8SToby Isaac . imageDim - the dimension of the image of the transformation 370*907761f8SToby Isaac . origin - a point in the original space 371*907761f8SToby Isaac . originImage - the image of the origin under the transformation 372*907761f8SToby Isaac . J - the Jacobian of the image: an [imageDim x dim] matrix in row major order 373*907761f8SToby Isaac - formIndex - transform the quadrature weights as k-forms of this index (if the number of components is a multiple of (dim choose formIndex), it is assumed that they represent multiple k-forms) [see PetscDTAltVPullback() for interpretation of formIndex] 374*907761f8SToby Isaac 375*907761f8SToby Isaac Output Arguments: 376*907761f8SToby Isaac . Jinvstarq - a quadrature rule where each point is the image of a point in the original quadrature rule, and where the k-form weights have been pulled-back by the pseudoinverse of J to the k-form weights in the image space. 377*907761f8SToby Isaac 378*907761f8SToby Isaac Note: the new quadrature rule will have a different number of components if spaces have different dimensions. For example, pushing a 2-form forward from a two dimensional space to a three dimensional space changes the number of components from 1 to 3. 379*907761f8SToby Isaac 380*907761f8SToby Isaac .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 381*907761f8SToby Isaac @*/ 382*907761f8SToby Isaac PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formIndex, PetscQuadrature *Jinvstarq) 383*907761f8SToby Isaac { 384*907761f8SToby Isaac PetscInt dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c; 385*907761f8SToby Isaac const PetscReal *points; 386*907761f8SToby Isaac const PetscReal *weights; 387*907761f8SToby Isaac PetscReal *imagePoints, *imageWeights; 388*907761f8SToby Isaac PetscReal *Jinv; 389*907761f8SToby Isaac PetscReal *Jinvstar; 390*907761f8SToby Isaac PetscErrorCode ierr; 391*907761f8SToby Isaac 392*907761f8SToby Isaac PetscFunctionBegin; 393*907761f8SToby Isaac PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 394*907761f8SToby Isaac if (imageDim < PetscAbsInt(formIndex)) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %D form in %D dimensions", PetscAbsInt(formIndex), imageDim); 395*907761f8SToby Isaac ierr = PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights);CHKERRQ(ierr); 396*907761f8SToby Isaac if (Nc % formSize) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Number of components %D is not a multiple of formSize %D\n", Nc, formSize); 397*907761f8SToby Isaac Ncopies = Nc / formSize; 398*907761f8SToby Isaac ierr = PetscDTBinomialInt(dim, PetscAbsInt(formIndex), &formSize);CHKERRQ(ierr); 399*907761f8SToby Isaac ierr = PetscDTBinomialInt(imageDim, PetscAbsInt(formIndex), &imageFormSize);CHKERRQ(ierr); 400*907761f8SToby Isaac imageNc = Ncopies * imageFormSize; 401*907761f8SToby Isaac ierr = PetscMalloc1(Npoints * imageDim, &imagePoints);CHKERRQ(ierr); 402*907761f8SToby Isaac ierr = PetscMalloc1(Npoints * imageNc, &imageWeights);CHKERRQ(ierr); 403*907761f8SToby Isaac ierr = PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar);CHKERRQ(ierr); 404*907761f8SToby Isaac ierr = PetscDTJacobianInverse_Internal(dim, imageDim, J, Jinv);CHKERRQ(ierr); 405*907761f8SToby Isaac ierr = PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formIndex, Jinvstar);CHKERRQ(ierr); 406*907761f8SToby Isaac for (pt = 0; pt < Npoints; pt++) { 407*907761f8SToby Isaac const PetscReal *point = &points[pt * dim]; 408*907761f8SToby Isaac PetscReal *imagePoint = &imagePoints[pt * imageDim]; 409*907761f8SToby Isaac 410*907761f8SToby Isaac for (i = 0; i < imageDim; i++) { 411*907761f8SToby Isaac PetscReal val = originImage[i]; 412*907761f8SToby Isaac 413*907761f8SToby Isaac for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]); 414*907761f8SToby Isaac imagePoint[i] = val; 415*907761f8SToby Isaac } 416*907761f8SToby Isaac for (c = 0; c < Ncopies; c++) { 417*907761f8SToby Isaac const PetscReal *form = &weights[pt * Nc + c * formSize]; 418*907761f8SToby Isaac PetscReal *imageForm = &imageWeights[pt * imageNc + c * imageFormSize]; 419*907761f8SToby Isaac 420*907761f8SToby Isaac for (i = 0; i < imageFormSize; i++) { 421*907761f8SToby Isaac PetscReal val = 0.; 422*907761f8SToby Isaac 423*907761f8SToby Isaac for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j]; 424*907761f8SToby Isaac imageForm[i] = val; 425*907761f8SToby Isaac } 426*907761f8SToby Isaac } 427*907761f8SToby Isaac } 428*907761f8SToby Isaac ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq);CHKERRQ(ierr); 429*907761f8SToby Isaac ierr = PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights);CHKERRQ(ierr); 430*907761f8SToby Isaac ierr = PetscFree2(Jinv, Jinvstar);CHKERRQ(ierr); 431*907761f8SToby Isaac PetscFunctionReturn(0); 432*907761f8SToby Isaac } 433*907761f8SToby Isaac 43440d8ff71SMatthew G. Knepley /*@C 43540d8ff71SMatthew G. Knepley PetscQuadratureSetData - Sets the data defining the quadrature 43640d8ff71SMatthew G. Knepley 43740d8ff71SMatthew G. Knepley Not collective 43840d8ff71SMatthew G. Knepley 43940d8ff71SMatthew G. Knepley Input Parameters: 44040d8ff71SMatthew G. Knepley + q - The PetscQuadrature object 44140d8ff71SMatthew G. Knepley . dim - The spatial dimension 442e2b35d93SBarry Smith . Nc - The number of components 44340d8ff71SMatthew G. Knepley . npoints - The number of quadrature points 44440d8ff71SMatthew G. Knepley . points - The coordinates of each quadrature point 44540d8ff71SMatthew G. Knepley - weights - The weight of each quadrature point 44640d8ff71SMatthew G. Knepley 447c99e0549SMatthew 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. 448f2fd9e53SMatthew G. Knepley 44940d8ff71SMatthew G. Knepley Level: intermediate 45040d8ff71SMatthew G. Knepley 45140d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 45240d8ff71SMatthew G. Knepley @*/ 453a6b92713SMatthew G. Knepley PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 45421454ff5SMatthew G. Knepley { 45521454ff5SMatthew G. Knepley PetscFunctionBegin; 45621454ff5SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 45721454ff5SMatthew G. Knepley if (dim >= 0) q->dim = dim; 458a6b92713SMatthew G. Knepley if (Nc >= 0) q->Nc = Nc; 45921454ff5SMatthew G. Knepley if (npoints >= 0) q->numPoints = npoints; 46021454ff5SMatthew G. Knepley if (points) { 46121454ff5SMatthew G. Knepley PetscValidPointer(points, 4); 46221454ff5SMatthew G. Knepley q->points = points; 46321454ff5SMatthew G. Knepley } 46421454ff5SMatthew G. Knepley if (weights) { 46521454ff5SMatthew G. Knepley PetscValidPointer(weights, 5); 46621454ff5SMatthew G. Knepley q->weights = weights; 46721454ff5SMatthew G. Knepley } 468f9fd7fdbSMatthew G. Knepley PetscFunctionReturn(0); 469f9fd7fdbSMatthew G. Knepley } 470f9fd7fdbSMatthew G. Knepley 471d9bac1caSLisandro Dalcin static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v) 472d9bac1caSLisandro Dalcin { 473d9bac1caSLisandro Dalcin PetscInt q, d, c; 474d9bac1caSLisandro Dalcin PetscViewerFormat format; 475d9bac1caSLisandro Dalcin PetscErrorCode ierr; 476d9bac1caSLisandro Dalcin 477d9bac1caSLisandro Dalcin PetscFunctionBegin; 478c74b4a09SMatthew 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);} 479c74b4a09SMatthew G. Knepley else {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);} 480d9bac1caSLisandro Dalcin ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 481d9bac1caSLisandro Dalcin if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0); 482d9bac1caSLisandro Dalcin for (q = 0; q < quad->numPoints; ++q) { 483c74b4a09SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr); 484d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr); 485d9bac1caSLisandro Dalcin for (d = 0; d < quad->dim; ++d) { 486d9bac1caSLisandro Dalcin if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 487d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 488d9bac1caSLisandro Dalcin } 489d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr); 490c74b4a09SMatthew G. Knepley if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);} 491d9bac1caSLisandro Dalcin for (c = 0; c < quad->Nc; ++c) { 492d9bac1caSLisandro Dalcin if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 493c74b4a09SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr); 494d9bac1caSLisandro Dalcin } 495d9bac1caSLisandro Dalcin if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);} 496d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr); 497d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr); 498d9bac1caSLisandro Dalcin } 499d9bac1caSLisandro Dalcin PetscFunctionReturn(0); 500d9bac1caSLisandro Dalcin } 501d9bac1caSLisandro Dalcin 50240d8ff71SMatthew G. Knepley /*@C 50340d8ff71SMatthew G. Knepley PetscQuadratureView - Views a PetscQuadrature object 50440d8ff71SMatthew G. Knepley 505d083f849SBarry Smith Collective on quad 50640d8ff71SMatthew G. Knepley 50740d8ff71SMatthew G. Knepley Input Parameters: 508d9bac1caSLisandro Dalcin + quad - The PetscQuadrature object 50940d8ff71SMatthew G. Knepley - viewer - The PetscViewer object 51040d8ff71SMatthew G. Knepley 51140d8ff71SMatthew G. Knepley Level: beginner 51240d8ff71SMatthew G. Knepley 51340d8ff71SMatthew G. Knepley .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 51440d8ff71SMatthew G. Knepley @*/ 515f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 516f9fd7fdbSMatthew G. Knepley { 517d9bac1caSLisandro Dalcin PetscBool iascii; 518f9fd7fdbSMatthew G. Knepley PetscErrorCode ierr; 519f9fd7fdbSMatthew G. Knepley 520f9fd7fdbSMatthew G. Knepley PetscFunctionBegin; 521d9bac1caSLisandro Dalcin PetscValidHeader(quad, 1); 522d9bac1caSLisandro Dalcin if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 523d9bac1caSLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);} 524d9bac1caSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 525d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 526d9bac1caSLisandro Dalcin if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);} 527d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 528bfa639d9SMatthew G. Knepley PetscFunctionReturn(0); 529bfa639d9SMatthew G. Knepley } 530bfa639d9SMatthew G. Knepley 53189710940SMatthew G. Knepley /*@C 53289710940SMatthew G. Knepley PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 53389710940SMatthew G. Knepley 53489710940SMatthew G. Knepley Not collective 53589710940SMatthew G. Knepley 53689710940SMatthew G. Knepley Input Parameter: 53789710940SMatthew G. Knepley + q - The original PetscQuadrature 53889710940SMatthew G. Knepley . numSubelements - The number of subelements the original element is divided into 53989710940SMatthew G. Knepley . v0 - An array of the initial points for each subelement 54089710940SMatthew G. Knepley - jac - An array of the Jacobian mappings from the reference to each subelement 54189710940SMatthew G. Knepley 54289710940SMatthew G. Knepley Output Parameters: 54389710940SMatthew G. Knepley . dim - The dimension 54489710940SMatthew G. Knepley 54589710940SMatthew G. Knepley Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 54689710940SMatthew G. Knepley 547f5f57ec0SBarry Smith Not available from Fortran 548f5f57ec0SBarry Smith 54989710940SMatthew G. Knepley Level: intermediate 55089710940SMatthew G. Knepley 55189710940SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 55289710940SMatthew G. Knepley @*/ 55389710940SMatthew G. Knepley PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 55489710940SMatthew G. Knepley { 55589710940SMatthew G. Knepley const PetscReal *points, *weights; 55689710940SMatthew G. Knepley PetscReal *pointsRef, *weightsRef; 557a6b92713SMatthew G. Knepley PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e; 55889710940SMatthew G. Knepley PetscErrorCode ierr; 55989710940SMatthew G. Knepley 56089710940SMatthew G. Knepley PetscFunctionBegin; 56189710940SMatthew G. Knepley PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 56289710940SMatthew G. Knepley PetscValidPointer(v0, 3); 56389710940SMatthew G. Knepley PetscValidPointer(jac, 4); 56489710940SMatthew G. Knepley PetscValidPointer(qref, 5); 56589710940SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 56689710940SMatthew G. Knepley ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 567a6b92713SMatthew G. Knepley ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 56889710940SMatthew G. Knepley npointsRef = npoints*numSubelements; 56989710940SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 570a6b92713SMatthew G. Knepley ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr); 57189710940SMatthew G. Knepley for (c = 0; c < numSubelements; ++c) { 57289710940SMatthew G. Knepley for (p = 0; p < npoints; ++p) { 57389710940SMatthew G. Knepley for (d = 0; d < dim; ++d) { 57489710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 57589710940SMatthew G. Knepley for (e = 0; e < dim; ++e) { 57689710940SMatthew G. Knepley pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 57789710940SMatthew G. Knepley } 57889710940SMatthew G. Knepley } 57989710940SMatthew G. Knepley /* Could also use detJ here */ 580a6b92713SMatthew G. Knepley for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements; 58189710940SMatthew G. Knepley } 58289710940SMatthew G. Knepley } 58389710940SMatthew G. Knepley ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 584a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 58589710940SMatthew G. Knepley PetscFunctionReturn(0); 58689710940SMatthew G. Knepley } 58789710940SMatthew G. Knepley 58837045ce4SJed Brown /*@ 58937045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 59037045ce4SJed Brown 59137045ce4SJed Brown Not Collective 59237045ce4SJed Brown 59337045ce4SJed Brown Input Arguments: 59437045ce4SJed Brown + npoints - number of spatial points to evaluate at 59537045ce4SJed Brown . points - array of locations to evaluate at 59637045ce4SJed Brown . ndegree - number of basis degrees to evaluate 59737045ce4SJed Brown - degrees - sorted array of degrees to evaluate 59837045ce4SJed Brown 59937045ce4SJed Brown Output Arguments: 6000298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 6010298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 6020298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 60337045ce4SJed Brown 60437045ce4SJed Brown Level: intermediate 60537045ce4SJed Brown 60637045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 60737045ce4SJed Brown @*/ 60837045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 60937045ce4SJed Brown { 61037045ce4SJed Brown PetscInt i,maxdegree; 61137045ce4SJed Brown 61237045ce4SJed Brown PetscFunctionBegin; 61337045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 61437045ce4SJed Brown maxdegree = degrees[ndegree-1]; 61537045ce4SJed Brown for (i=0; i<npoints; i++) { 61637045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 61737045ce4SJed Brown PetscInt j,k; 61837045ce4SJed Brown x = points[i]; 61937045ce4SJed Brown pm2 = 0; 62037045ce4SJed Brown pm1 = 1; 62137045ce4SJed Brown pd2 = 0; 62237045ce4SJed Brown pd1 = 0; 62337045ce4SJed Brown pdd2 = 0; 62437045ce4SJed Brown pdd1 = 0; 62537045ce4SJed Brown k = 0; 62637045ce4SJed Brown if (degrees[k] == 0) { 62737045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 62837045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 62937045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 63037045ce4SJed Brown k++; 63137045ce4SJed Brown } 63237045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 63337045ce4SJed Brown PetscReal p,d,dd; 63437045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 63537045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 63637045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 63737045ce4SJed Brown pm2 = pm1; 63837045ce4SJed Brown pm1 = p; 63937045ce4SJed Brown pd2 = pd1; 64037045ce4SJed Brown pd1 = d; 64137045ce4SJed Brown pdd2 = pdd1; 64237045ce4SJed Brown pdd1 = dd; 64337045ce4SJed Brown if (degrees[k] == j) { 64437045ce4SJed Brown if (B) B[i*ndegree+k] = p; 64537045ce4SJed Brown if (D) D[i*ndegree+k] = d; 64637045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 64737045ce4SJed Brown } 64837045ce4SJed Brown } 64937045ce4SJed Brown } 65037045ce4SJed Brown PetscFunctionReturn(0); 65137045ce4SJed Brown } 65237045ce4SJed Brown 65337045ce4SJed Brown /*@ 65437045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 65537045ce4SJed Brown 65637045ce4SJed Brown Not Collective 65737045ce4SJed Brown 65837045ce4SJed Brown Input Arguments: 65937045ce4SJed Brown + npoints - number of points 66037045ce4SJed Brown . a - left end of interval (often-1) 66137045ce4SJed Brown - b - right end of interval (often +1) 66237045ce4SJed Brown 66337045ce4SJed Brown Output Arguments: 66437045ce4SJed Brown + x - quadrature points 66537045ce4SJed Brown - w - quadrature weights 66637045ce4SJed Brown 66737045ce4SJed Brown Level: intermediate 66837045ce4SJed Brown 66937045ce4SJed Brown References: 67096a0c994SBarry Smith . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 67137045ce4SJed Brown 67237045ce4SJed Brown .seealso: PetscDTLegendreEval() 67337045ce4SJed Brown @*/ 67437045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 67537045ce4SJed Brown { 67637045ce4SJed Brown PetscErrorCode ierr; 67737045ce4SJed Brown PetscInt i; 67837045ce4SJed Brown PetscReal *work; 67937045ce4SJed Brown PetscScalar *Z; 68037045ce4SJed Brown PetscBLASInt N,LDZ,info; 68137045ce4SJed Brown 68237045ce4SJed Brown PetscFunctionBegin; 6830bfcf5a5SMatthew G. Knepley ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 68437045ce4SJed Brown /* Set up the Golub-Welsch system */ 68537045ce4SJed Brown for (i=0; i<npoints; i++) { 68637045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 68737045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 68837045ce4SJed Brown } 689dcca6d9dSJed Brown ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 690c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 69137045ce4SJed Brown LDZ = N; 69237045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 6938b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 69437045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 6951c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 69637045ce4SJed Brown 69737045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 69837045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 69937045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 70019a57d60SBarry Smith if (x[i] == -0.0) x[i] = 0.0; 70137045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 7020d644c17SKarl Rupp 70388393a60SJed Brown w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 70437045ce4SJed Brown } 70537045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 70637045ce4SJed Brown PetscFunctionReturn(0); 70737045ce4SJed Brown } 708194825f6SJed Brown 7098272889dSSatish Balay static void qAndLEvaluation(PetscInt n, PetscReal x, PetscReal *q, PetscReal *qp, PetscReal *Ln) 7108272889dSSatish Balay /* 7118272889dSSatish Balay Compute the polynomial q(x) = L_{N+1}(x) - L_{n-1}(x) and its derivative in 7128272889dSSatish Balay addition to L_N(x) as these are needed for computing the GLL points via Newton's method. 7138272889dSSatish Balay Reference: "Implementing Spectral Methods for Partial Differential Equations: Algorithms 7148272889dSSatish Balay for Scientists and Engineers" by David A. Kopriva. 7158272889dSSatish Balay */ 7168272889dSSatish Balay { 7178272889dSSatish Balay PetscInt k; 7188272889dSSatish Balay 7198272889dSSatish Balay PetscReal Lnp; 7208272889dSSatish Balay PetscReal Lnp1, Lnp1p; 7218272889dSSatish Balay PetscReal Lnm1, Lnm1p; 7228272889dSSatish Balay PetscReal Lnm2, Lnm2p; 7238272889dSSatish Balay 7248272889dSSatish Balay Lnm1 = 1.0; 7258272889dSSatish Balay *Ln = x; 7268272889dSSatish Balay Lnm1p = 0.0; 7278272889dSSatish Balay Lnp = 1.0; 7288272889dSSatish Balay 7298272889dSSatish Balay for (k=2; k<=n; ++k) { 7308272889dSSatish Balay Lnm2 = Lnm1; 7318272889dSSatish Balay Lnm1 = *Ln; 7328272889dSSatish Balay Lnm2p = Lnm1p; 7338272889dSSatish Balay Lnm1p = Lnp; 7348272889dSSatish Balay *Ln = (2.*((PetscReal)k)-1.)/(1.0*((PetscReal)k))*x*Lnm1 - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm2; 7358272889dSSatish Balay Lnp = Lnm2p + (2.0*((PetscReal)k)-1.)*Lnm1; 7368272889dSSatish Balay } 7378272889dSSatish Balay k = n+1; 7388272889dSSatish Balay Lnp1 = (2.*((PetscReal)k)-1.)/(((PetscReal)k))*x*(*Ln) - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm1; 7398272889dSSatish Balay Lnp1p = Lnm1p + (2.0*((PetscReal)k)-1.)*(*Ln); 7408272889dSSatish Balay *q = Lnp1 - Lnm1; 7418272889dSSatish Balay *qp = Lnp1p - Lnm1p; 7428272889dSSatish Balay } 7438272889dSSatish Balay 7448272889dSSatish Balay /*@C 7458272889dSSatish Balay PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre 7468272889dSSatish Balay nodes of a given size on the domain [-1,1] 7478272889dSSatish Balay 7488272889dSSatish Balay Not Collective 7498272889dSSatish Balay 7508272889dSSatish Balay Input Parameter: 7518272889dSSatish Balay + n - number of grid nodes 752f2e8fe4dShannah_mairs - type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 7538272889dSSatish Balay 7548272889dSSatish Balay Output Arguments: 7558272889dSSatish Balay + x - quadrature points 7568272889dSSatish Balay - w - quadrature weights 7578272889dSSatish Balay 7588272889dSSatish Balay Notes: 7598272889dSSatish Balay For n > 30 the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not 7608272889dSSatish Balay close enough to the desired solution 7618272889dSSatish Balay 7628272889dSSatish Balay These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes 7638272889dSSatish Balay 764a8d69d7bSBarry 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 7658272889dSSatish Balay 7668272889dSSatish Balay Level: intermediate 7678272889dSSatish Balay 7688272889dSSatish Balay .seealso: PetscDTGaussQuadrature() 7698272889dSSatish Balay 7708272889dSSatish Balay @*/ 771916e780bShannah_mairs PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w) 7728272889dSSatish Balay { 7738272889dSSatish Balay PetscErrorCode ierr; 7748272889dSSatish Balay 7758272889dSSatish Balay PetscFunctionBegin; 7768272889dSSatish Balay if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element"); 7778272889dSSatish Balay 778f2e8fe4dShannah_mairs if (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA) { 7798272889dSSatish Balay PetscReal *M,si; 7808272889dSSatish Balay PetscBLASInt bn,lierr; 7818272889dSSatish Balay PetscReal x0,z0,z1,z2; 7828272889dSSatish Balay PetscInt i,p = npoints - 1,nn; 7838272889dSSatish Balay 7848272889dSSatish Balay x[0] =-1.0; 7858272889dSSatish Balay x[npoints-1] = 1.0; 7868272889dSSatish Balay if (npoints-2 > 0){ 7878272889dSSatish Balay ierr = PetscMalloc1(npoints-1,&M);CHKERRQ(ierr); 7888272889dSSatish Balay for (i=0; i<npoints-2; i++) { 7898272889dSSatish Balay si = ((PetscReal)i)+1.0; 7908272889dSSatish Balay M[i]=0.5*PetscSqrtReal(si*(si+2.0)/((si+0.5)*(si+1.5))); 7918272889dSSatish Balay } 7928272889dSSatish Balay ierr = PetscBLASIntCast(npoints-2,&bn);CHKERRQ(ierr); 793580bdb30SBarry Smith ierr = PetscArrayzero(&x[1],bn);CHKERRQ(ierr); 7948272889dSSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 7958272889dSSatish Balay x0=0; 7968272889dSSatish Balay PetscStackCallBLAS("LAPACKsteqr",LAPACKREALsteqr_("N",&bn,&x[1],M,&x0,&bn,M,&lierr)); 7978272889dSSatish Balay if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in STERF Lapack routine %d",(int)lierr); 7988272889dSSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 7998272889dSSatish Balay ierr = PetscFree(M);CHKERRQ(ierr); 8008272889dSSatish Balay } 8018272889dSSatish Balay if ((npoints-1)%2==0) { 8028272889dSSatish Balay x[(npoints-1)/2] = 0.0; /* hard wire to exactly 0.0 since linear algebra produces nonzero */ 8038272889dSSatish Balay } 8048272889dSSatish Balay 8058272889dSSatish Balay w[0] = w[p] = 2.0/(((PetscReal)(p))*(((PetscReal)p)+1.0)); 8068272889dSSatish Balay z2 = -1.; /* Dummy value to avoid -Wmaybe-initialized */ 8078272889dSSatish Balay for (i=1; i<p; i++) { 8088272889dSSatish Balay x0 = x[i]; 8098272889dSSatish Balay z0 = 1.0; 8108272889dSSatish Balay z1 = x0; 8118272889dSSatish Balay for (nn=1; nn<p; nn++) { 8128272889dSSatish Balay z2 = x0*z1*(2.0*((PetscReal)nn)+1.0)/(((PetscReal)nn)+1.0)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.0)); 8138272889dSSatish Balay z0 = z1; 8148272889dSSatish Balay z1 = z2; 8158272889dSSatish Balay } 8168272889dSSatish Balay w[i]=2.0/(((PetscReal)p)*(((PetscReal)p)+1.0)*z2*z2); 8178272889dSSatish Balay } 8188272889dSSatish Balay } else { 8198272889dSSatish Balay PetscInt j,m; 8208272889dSSatish Balay PetscReal z1,z,q,qp,Ln; 8218272889dSSatish Balay PetscReal *pt; 8228272889dSSatish Balay ierr = PetscMalloc1(npoints,&pt);CHKERRQ(ierr); 8238272889dSSatish Balay 824d410ae54Shannah_mairs if (npoints > 30) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON produces incorrect answers for n > 30"); 8258272889dSSatish Balay x[0] = -1.0; 8268272889dSSatish Balay x[npoints-1] = 1.0; 827feb237baSPierre Jolivet w[0] = w[npoints-1] = 2./(((PetscReal)npoints)*(((PetscReal)npoints)-1.0)); 8288272889dSSatish Balay m = (npoints-1)/2; /* The roots are symmetric, so we only find half of them. */ 8298272889dSSatish Balay for (j=1; j<=m; j++) { /* Loop over the desired roots. */ 8308272889dSSatish 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)); 8318272889dSSatish Balay /* Starting with the above approximation to the ith root, we enter */ 8328272889dSSatish Balay /* the main loop of refinement by Newton's method. */ 8338272889dSSatish Balay do { 8348272889dSSatish Balay qAndLEvaluation(npoints-1,z,&q,&qp,&Ln); 8358272889dSSatish Balay z1 = z; 8368272889dSSatish Balay z = z1-q/qp; /* Newton's method. */ 8378272889dSSatish Balay } while (PetscAbs(z-z1) > 10.*PETSC_MACHINE_EPSILON); 8388272889dSSatish Balay qAndLEvaluation(npoints-1,z,&q,&qp,&Ln); 8398272889dSSatish Balay 8408272889dSSatish Balay x[j] = z; 8418272889dSSatish Balay x[npoints-1-j] = -z; /* and put in its symmetric counterpart. */ 8428272889dSSatish Balay w[j] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln); /* Compute the weight */ 8438272889dSSatish Balay w[npoints-1-j] = w[j]; /* and its symmetric counterpart. */ 8448272889dSSatish Balay pt[j]=qp; 8458272889dSSatish Balay } 8468272889dSSatish Balay 8478272889dSSatish Balay if ((npoints-1)%2==0) { 8488272889dSSatish Balay qAndLEvaluation(npoints-1,0.0,&q,&qp,&Ln); 8498272889dSSatish Balay x[(npoints-1)/2] = 0.0; 8508272889dSSatish Balay w[(npoints-1)/2] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln); 8518272889dSSatish Balay } 8528272889dSSatish Balay ierr = PetscFree(pt);CHKERRQ(ierr); 8538272889dSSatish Balay } 8548272889dSSatish Balay PetscFunctionReturn(0); 8558272889dSSatish Balay } 8568272889dSSatish Balay 857744bafbcSMatthew G. Knepley /*@ 858744bafbcSMatthew G. Knepley PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 859744bafbcSMatthew G. Knepley 860744bafbcSMatthew G. Knepley Not Collective 861744bafbcSMatthew G. Knepley 862744bafbcSMatthew G. Knepley Input Arguments: 863744bafbcSMatthew G. Knepley + dim - The spatial dimension 864a6b92713SMatthew G. Knepley . Nc - The number of components 865744bafbcSMatthew G. Knepley . npoints - number of points in one dimension 866744bafbcSMatthew G. Knepley . a - left end of interval (often-1) 867744bafbcSMatthew G. Knepley - b - right end of interval (often +1) 868744bafbcSMatthew G. Knepley 869744bafbcSMatthew G. Knepley Output Argument: 870744bafbcSMatthew G. Knepley . q - A PetscQuadrature object 871744bafbcSMatthew G. Knepley 872744bafbcSMatthew G. Knepley Level: intermediate 873744bafbcSMatthew G. Knepley 874744bafbcSMatthew G. Knepley .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 875744bafbcSMatthew G. Knepley @*/ 876a6b92713SMatthew G. Knepley PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 877744bafbcSMatthew G. Knepley { 878a6b92713SMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c; 879744bafbcSMatthew G. Knepley PetscReal *x, *w, *xw, *ww; 880744bafbcSMatthew G. Knepley PetscErrorCode ierr; 881744bafbcSMatthew G. Knepley 882744bafbcSMatthew G. Knepley PetscFunctionBegin; 883744bafbcSMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 884a6b92713SMatthew G. Knepley ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr); 885744bafbcSMatthew G. Knepley /* Set up the Golub-Welsch system */ 886744bafbcSMatthew G. Knepley switch (dim) { 887744bafbcSMatthew G. Knepley case 0: 888744bafbcSMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 889744bafbcSMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 890744bafbcSMatthew G. Knepley ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 891a6b92713SMatthew G. Knepley ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 892744bafbcSMatthew G. Knepley x[0] = 0.0; 893a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[c] = 1.0; 894744bafbcSMatthew G. Knepley break; 895744bafbcSMatthew G. Knepley case 1: 896a6b92713SMatthew G. Knepley ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr); 897a6b92713SMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr); 898a6b92713SMatthew G. Knepley for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i]; 899a6b92713SMatthew G. Knepley ierr = PetscFree(ww);CHKERRQ(ierr); 900744bafbcSMatthew G. Knepley break; 901744bafbcSMatthew G. Knepley case 2: 902744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 903744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 904744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 905744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 906744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+0] = xw[i]; 907744bafbcSMatthew G. Knepley x[(i*npoints+j)*dim+1] = xw[j]; 908a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j]; 909744bafbcSMatthew G. Knepley } 910744bafbcSMatthew G. Knepley } 911744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 912744bafbcSMatthew G. Knepley break; 913744bafbcSMatthew G. Knepley case 3: 914744bafbcSMatthew G. Knepley ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 915744bafbcSMatthew G. Knepley ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 916744bafbcSMatthew G. Knepley for (i = 0; i < npoints; ++i) { 917744bafbcSMatthew G. Knepley for (j = 0; j < npoints; ++j) { 918744bafbcSMatthew G. Knepley for (k = 0; k < npoints; ++k) { 919744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 920744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 921744bafbcSMatthew G. Knepley x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 922a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k]; 923744bafbcSMatthew G. Knepley } 924744bafbcSMatthew G. Knepley } 925744bafbcSMatthew G. Knepley } 926744bafbcSMatthew G. Knepley ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 927744bafbcSMatthew G. Knepley break; 928744bafbcSMatthew G. Knepley default: 929744bafbcSMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 930744bafbcSMatthew G. Knepley } 931744bafbcSMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 9322f5fb066SToby Isaac ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 933a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 934d9bac1caSLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr); 935744bafbcSMatthew G. Knepley PetscFunctionReturn(0); 936744bafbcSMatthew G. Knepley } 937744bafbcSMatthew G. Knepley 938494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 939494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 940494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 941494e7359SMatthew G. Knepley { 942494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 943494e7359SMatthew G. Knepley PetscInt k; 944494e7359SMatthew G. Knepley 945494e7359SMatthew G. Knepley PetscFunctionBegin; 946494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 947494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 948494e7359SMatthew G. Knepley apb = a + b; 949494e7359SMatthew G. Knepley pn2 = 1.0; 950494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 951494e7359SMatthew G. Knepley *P = 0.0; 952494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 953494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 954494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 955494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 956494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 957494e7359SMatthew G. Knepley 958494e7359SMatthew G. Knepley a2 = a2 / a1; 959494e7359SMatthew G. Knepley a3 = a3 / a1; 960494e7359SMatthew G. Knepley a4 = a4 / a1; 961494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 962494e7359SMatthew G. Knepley pn2 = pn1; 963494e7359SMatthew G. Knepley pn1 = *P; 964494e7359SMatthew G. Knepley } 965494e7359SMatthew G. Knepley PetscFunctionReturn(0); 966494e7359SMatthew G. Knepley } 967494e7359SMatthew G. Knepley 968494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 969494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 970494e7359SMatthew G. Knepley { 971494e7359SMatthew G. Knepley PetscReal nP; 972494e7359SMatthew G. Knepley PetscErrorCode ierr; 973494e7359SMatthew G. Knepley 974494e7359SMatthew G. Knepley PetscFunctionBegin; 975494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 976494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 977494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 978494e7359SMatthew G. Knepley PetscFunctionReturn(0); 979494e7359SMatthew G. Knepley } 980494e7359SMatthew G. Knepley 981494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 982494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 983494e7359SMatthew G. Knepley { 984494e7359SMatthew G. Knepley PetscFunctionBegin; 985494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 986494e7359SMatthew G. Knepley *eta = y; 987494e7359SMatthew G. Knepley PetscFunctionReturn(0); 988494e7359SMatthew G. Knepley } 989494e7359SMatthew G. Knepley 990494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 991494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 992494e7359SMatthew G. Knepley { 993494e7359SMatthew G. Knepley PetscFunctionBegin; 994494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 995494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 996494e7359SMatthew G. Knepley *zeta = z; 997494e7359SMatthew G. Knepley PetscFunctionReturn(0); 998494e7359SMatthew G. Knepley } 999494e7359SMatthew G. Knepley 1000494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 1001494e7359SMatthew G. Knepley { 1002494e7359SMatthew G. Knepley PetscInt maxIter = 100; 1003494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 1004a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 1005494e7359SMatthew G. Knepley PetscInt k; 1006494e7359SMatthew G. Knepley PetscErrorCode ierr; 1007494e7359SMatthew G. Knepley 1008494e7359SMatthew G. Knepley PetscFunctionBegin; 1009a8291ba1SSatish Balay 10108b49ba18SBarry Smith a1 = PetscPowReal(2.0, a+b+1); 1011a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 10120646a658SBarry Smith a2 = PetscTGamma(a + npoints + 1); 10130646a658SBarry Smith a3 = PetscTGamma(b + npoints + 1); 10140646a658SBarry Smith a4 = PetscTGamma(a + b + npoints + 1); 1015a8291ba1SSatish Balay #else 101629bcbfd0SToby Isaac { 1017d24bbb91SToby Isaac PetscInt ia, ib; 101829bcbfd0SToby Isaac 1019d24bbb91SToby Isaac ia = (PetscInt) a; 1020d24bbb91SToby Isaac ib = (PetscInt) b; 1021d24bbb91SToby 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 */ 1022fad4db65SToby Isaac ierr = PetscDTFactorial(ia + npoints, &a2);CHKERRQ(ierr); 1023fad4db65SToby Isaac ierr = PetscDTFactorial(ib + npoints, &a3);CHKERRQ(ierr); 1024fad4db65SToby Isaac ierr = PetscDTFactorial(ia + ib + npoints, &a4);CHKERRQ(ierr); 102529bcbfd0SToby Isaac } else { 1026a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 102729bcbfd0SToby Isaac } 102829bcbfd0SToby Isaac } 1029a8291ba1SSatish Balay #endif 1030a8291ba1SSatish Balay 1031fad4db65SToby Isaac ierr = PetscDTFactorial(npoints, &a5);CHKERRQ(ierr); 1032494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 1033494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 1034494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 1035494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 10368b49ba18SBarry Smith PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 1037494e7359SMatthew G. Knepley PetscInt j; 1038494e7359SMatthew G. Knepley 1039494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 1040494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 1041494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 1042494e7359SMatthew G. Knepley PetscInt i; 1043494e7359SMatthew G. Knepley 1044494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 1045494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 1046494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 1047494e7359SMatthew G. Knepley delta = f / (fp - f * s); 1048494e7359SMatthew G. Knepley r = r - delta; 104977b4d14cSPeter Brune if (PetscAbsReal(delta) < eps) break; 1050494e7359SMatthew G. Knepley } 1051494e7359SMatthew G. Knepley x[k] = r; 1052494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 1053494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 1054494e7359SMatthew G. Knepley } 1055494e7359SMatthew G. Knepley PetscFunctionReturn(0); 1056494e7359SMatthew G. Knepley } 1057494e7359SMatthew G. Knepley 1058f5f57ec0SBarry Smith /*@ 1059494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 1060494e7359SMatthew G. Knepley 1061494e7359SMatthew G. Knepley Not Collective 1062494e7359SMatthew G. Knepley 1063494e7359SMatthew G. Knepley Input Arguments: 1064494e7359SMatthew G. Knepley + dim - The simplex dimension 1065a6b92713SMatthew G. Knepley . Nc - The number of components 1066dcce0ee2SMatthew G. Knepley . npoints - The number of points in one dimension 1067494e7359SMatthew G. Knepley . a - left end of interval (often-1) 1068494e7359SMatthew G. Knepley - b - right end of interval (often +1) 1069494e7359SMatthew G. Knepley 1070744bafbcSMatthew G. Knepley Output Argument: 1071552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object 1072494e7359SMatthew G. Knepley 1073494e7359SMatthew G. Knepley Level: intermediate 1074494e7359SMatthew G. Knepley 1075494e7359SMatthew G. Knepley References: 107696a0c994SBarry Smith . 1. - Karniadakis and Sherwin. FIAT 1077494e7359SMatthew G. Knepley 1078744bafbcSMatthew G. Knepley .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 1079494e7359SMatthew G. Knepley @*/ 1080dcce0ee2SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 1081494e7359SMatthew G. Knepley { 1082dcce0ee2SMatthew G. Knepley PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints; 1083494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 1084a6b92713SMatthew G. Knepley PetscInt i, j, k, c; 1085494e7359SMatthew G. Knepley PetscErrorCode ierr; 1086494e7359SMatthew G. Knepley 1087494e7359SMatthew G. Knepley PetscFunctionBegin; 1088494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 1089dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr); 1090dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr); 1091494e7359SMatthew G. Knepley switch (dim) { 1092707aa5c5SMatthew G. Knepley case 0: 1093707aa5c5SMatthew G. Knepley ierr = PetscFree(x);CHKERRQ(ierr); 1094707aa5c5SMatthew G. Knepley ierr = PetscFree(w);CHKERRQ(ierr); 1095785e854fSJed Brown ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 1096a6b92713SMatthew G. Knepley ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 1097707aa5c5SMatthew G. Knepley x[0] = 0.0; 1098a6b92713SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[c] = 1.0; 1099707aa5c5SMatthew G. Knepley break; 1100494e7359SMatthew G. Knepley case 1: 1101dcce0ee2SMatthew G. Knepley ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr); 1102dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr); 1103dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i]; 1104a6b92713SMatthew G. Knepley ierr = PetscFree(wx);CHKERRQ(ierr); 1105494e7359SMatthew G. Knepley break; 1106494e7359SMatthew G. Knepley case 2: 1107dcce0ee2SMatthew G. Knepley ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr); 1108dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 1109dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 1110dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 1111dcce0ee2SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 1112dcce0ee2SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 1113dcce0ee2SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j]; 1114494e7359SMatthew G. Knepley } 1115494e7359SMatthew G. Knepley } 1116494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 1117494e7359SMatthew G. Knepley break; 1118494e7359SMatthew G. Knepley case 3: 1119dcce0ee2SMatthew G. Knepley ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr); 1120dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 1121dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 1122dcce0ee2SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 1123dcce0ee2SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 1124dcce0ee2SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 1125dcce0ee2SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 1126dcce0ee2SMatthew 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); 1127dcce0ee2SMatthew G. Knepley for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k]; 1128494e7359SMatthew G. Knepley } 1129494e7359SMatthew G. Knepley } 1130494e7359SMatthew G. Knepley } 1131494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 1132494e7359SMatthew G. Knepley break; 1133494e7359SMatthew G. Knepley default: 1134494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 1135494e7359SMatthew G. Knepley } 113621454ff5SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 11372f5fb066SToby Isaac ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 1138dcce0ee2SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 1139d9bac1caSLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr); 1140494e7359SMatthew G. Knepley PetscFunctionReturn(0); 1141494e7359SMatthew G. Knepley } 1142494e7359SMatthew G. Knepley 1143f5f57ec0SBarry Smith /*@ 1144b3c0f97bSTom Klotz PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 1145b3c0f97bSTom Klotz 1146b3c0f97bSTom Klotz Not Collective 1147b3c0f97bSTom Klotz 1148b3c0f97bSTom Klotz Input Arguments: 1149b3c0f97bSTom Klotz + dim - The cell dimension 1150b3c0f97bSTom Klotz . level - The number of points in one dimension, 2^l 1151b3c0f97bSTom Klotz . a - left end of interval (often-1) 1152b3c0f97bSTom Klotz - b - right end of interval (often +1) 1153b3c0f97bSTom Klotz 1154b3c0f97bSTom Klotz Output Argument: 1155b3c0f97bSTom Klotz . q - A PetscQuadrature object 1156b3c0f97bSTom Klotz 1157b3c0f97bSTom Klotz Level: intermediate 1158b3c0f97bSTom Klotz 1159b3c0f97bSTom Klotz .seealso: PetscDTGaussTensorQuadrature() 1160b3c0f97bSTom Klotz @*/ 1161b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 1162b3c0f97bSTom Klotz { 1163b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 1164b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1165b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1166b3c0f97bSTom Klotz const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 1167d84b4d08SMatthew G. Knepley PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 1168b3c0f97bSTom Klotz PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 1169b3c0f97bSTom Klotz PetscReal *x, *w; 1170b3c0f97bSTom Klotz PetscInt K, k, npoints; 1171b3c0f97bSTom Klotz PetscErrorCode ierr; 1172b3c0f97bSTom Klotz 1173b3c0f97bSTom Klotz PetscFunctionBegin; 1174b3c0f97bSTom Klotz if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 1175b3c0f97bSTom Klotz if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 1176b3c0f97bSTom Klotz /* Find K such that the weights are < 32 digits of precision */ 1177b3c0f97bSTom Klotz for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 11789add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 1179b3c0f97bSTom Klotz } 1180b3c0f97bSTom Klotz ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1181b3c0f97bSTom Klotz ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 1182b3c0f97bSTom Klotz npoints = 2*K-1; 1183b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 1184b3c0f97bSTom Klotz ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 1185b3c0f97bSTom Klotz /* Center term */ 1186b3c0f97bSTom Klotz x[0] = beta; 1187b3c0f97bSTom Klotz w[0] = 0.5*alpha*PETSC_PI; 1188b3c0f97bSTom Klotz for (k = 1; k < K; ++k) { 11899add2064SThomas Klotz wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 11901118d4bcSLisandro Dalcin xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h)); 1191b3c0f97bSTom Klotz x[2*k-1] = -alpha*xk+beta; 1192b3c0f97bSTom Klotz w[2*k-1] = wk; 1193b3c0f97bSTom Klotz x[2*k+0] = alpha*xk+beta; 1194b3c0f97bSTom Klotz w[2*k+0] = wk; 1195b3c0f97bSTom Klotz } 1196a6b92713SMatthew G. Knepley ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr); 1197b3c0f97bSTom Klotz PetscFunctionReturn(0); 1198b3c0f97bSTom Klotz } 1199b3c0f97bSTom Klotz 1200b3c0f97bSTom Klotz PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1201b3c0f97bSTom Klotz { 1202b3c0f97bSTom Klotz const PetscInt p = 16; /* Digits of precision in the evaluation */ 1203b3c0f97bSTom Klotz const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1204b3c0f97bSTom Klotz const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1205b3c0f97bSTom Klotz PetscReal h = 1.0; /* Step size, length between x_k */ 1206b3c0f97bSTom Klotz PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 1207b3c0f97bSTom Klotz PetscReal osum = 0.0; /* Integral on last level */ 1208b3c0f97bSTom Klotz PetscReal psum = 0.0; /* Integral on the level before the last level */ 1209b3c0f97bSTom Klotz PetscReal sum; /* Integral on current level */ 1210446c295cSMatthew G. Knepley PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 1211b3c0f97bSTom Klotz PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 1212b3c0f97bSTom Klotz PetscReal wk; /* Quadrature weight at x_k */ 1213b3c0f97bSTom Klotz PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 1214b3c0f97bSTom Klotz PetscInt d; /* Digits of precision in the integral */ 1215b3c0f97bSTom Klotz 1216b3c0f97bSTom Klotz PetscFunctionBegin; 1217b3c0f97bSTom Klotz if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 1218b3c0f97bSTom Klotz /* Center term */ 1219b3c0f97bSTom Klotz func(beta, &lval); 1220b3c0f97bSTom Klotz sum = 0.5*alpha*PETSC_PI*lval; 1221b3c0f97bSTom Klotz /* */ 1222b3c0f97bSTom Klotz do { 1223b3c0f97bSTom Klotz PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 1224b3c0f97bSTom Klotz PetscInt k = 1; 1225b3c0f97bSTom Klotz 1226b3c0f97bSTom Klotz ++l; 1227b3c0f97bSTom Klotz /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 1228b3c0f97bSTom Klotz /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 1229b3c0f97bSTom Klotz psum = osum; 1230b3c0f97bSTom Klotz osum = sum; 1231b3c0f97bSTom Klotz h *= 0.5; 1232b3c0f97bSTom Klotz sum *= 0.5; 1233b3c0f97bSTom Klotz do { 12349add2064SThomas Klotz wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1235446c295cSMatthew G. Knepley yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1236446c295cSMatthew G. Knepley lx = -alpha*(1.0 - yk)+beta; 1237446c295cSMatthew G. Knepley rx = alpha*(1.0 - yk)+beta; 1238b3c0f97bSTom Klotz func(lx, &lval); 1239b3c0f97bSTom Klotz func(rx, &rval); 1240b3c0f97bSTom Klotz lterm = alpha*wk*lval; 1241b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 1242b3c0f97bSTom Klotz sum += lterm; 1243b3c0f97bSTom Klotz rterm = alpha*wk*rval; 1244b3c0f97bSTom Klotz maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 1245b3c0f97bSTom Klotz sum += rterm; 1246b3c0f97bSTom Klotz ++k; 1247b3c0f97bSTom Klotz /* Only need to evaluate every other point on refined levels */ 1248b3c0f97bSTom Klotz if (l != 1) ++k; 12499add2064SThomas Klotz } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 1250b3c0f97bSTom Klotz 1251b3c0f97bSTom Klotz d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 1252b3c0f97bSTom Klotz d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 1253b3c0f97bSTom Klotz d3 = PetscLog10Real(maxTerm) - p; 125409d48545SBarry Smith if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 125509d48545SBarry Smith else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 1256b3c0f97bSTom Klotz d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 12579add2064SThomas Klotz } while (d < digits && l < 12); 1258b3c0f97bSTom Klotz *sol = sum; 1259e510cb1fSThomas Klotz 1260b3c0f97bSTom Klotz PetscFunctionReturn(0); 1261b3c0f97bSTom Klotz } 1262b3c0f97bSTom Klotz 1263497880caSRichard Tran Mills #if defined(PETSC_HAVE_MPFR) 126429f144ccSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 126529f144ccSMatthew G. Knepley { 1266e510cb1fSThomas Klotz const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 126729f144ccSMatthew G. Knepley PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 126829f144ccSMatthew G. Knepley mpfr_t alpha; /* Half-width of the integration interval */ 126929f144ccSMatthew G. Knepley mpfr_t beta; /* Center of the integration interval */ 127029f144ccSMatthew G. Knepley mpfr_t h; /* Step size, length between x_k */ 127129f144ccSMatthew G. Knepley mpfr_t osum; /* Integral on last level */ 127229f144ccSMatthew G. Knepley mpfr_t psum; /* Integral on the level before the last level */ 127329f144ccSMatthew G. Knepley mpfr_t sum; /* Integral on current level */ 127429f144ccSMatthew G. Knepley mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 127529f144ccSMatthew G. Knepley mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 127629f144ccSMatthew G. Knepley mpfr_t wk; /* Quadrature weight at x_k */ 127729f144ccSMatthew G. Knepley PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 127829f144ccSMatthew G. Knepley PetscInt d; /* Digits of precision in the integral */ 127929f144ccSMatthew G. Knepley mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 128029f144ccSMatthew G. Knepley 128129f144ccSMatthew G. Knepley PetscFunctionBegin; 128229f144ccSMatthew G. Knepley if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 128329f144ccSMatthew G. Knepley /* Create high precision storage */ 1284c9f744b5SMatthew 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); 128529f144ccSMatthew G. Knepley /* Initialization */ 128629f144ccSMatthew G. Knepley mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 128729f144ccSMatthew G. Knepley mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 128829f144ccSMatthew G. Knepley mpfr_set_d(osum, 0.0, MPFR_RNDN); 128929f144ccSMatthew G. Knepley mpfr_set_d(psum, 0.0, MPFR_RNDN); 129029f144ccSMatthew G. Knepley mpfr_set_d(h, 1.0, MPFR_RNDN); 129129f144ccSMatthew G. Knepley mpfr_const_pi(pi2, MPFR_RNDN); 129229f144ccSMatthew G. Knepley mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 129329f144ccSMatthew G. Knepley /* Center term */ 129429f144ccSMatthew G. Knepley func(0.5*(b+a), &lval); 129529f144ccSMatthew G. Knepley mpfr_set(sum, pi2, MPFR_RNDN); 129629f144ccSMatthew G. Knepley mpfr_mul(sum, sum, alpha, MPFR_RNDN); 129729f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 129829f144ccSMatthew G. Knepley /* */ 129929f144ccSMatthew G. Knepley do { 130029f144ccSMatthew G. Knepley PetscReal d1, d2, d3, d4; 130129f144ccSMatthew G. Knepley PetscInt k = 1; 130229f144ccSMatthew G. Knepley 130329f144ccSMatthew G. Knepley ++l; 130429f144ccSMatthew G. Knepley mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 130529f144ccSMatthew G. Knepley /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 130629f144ccSMatthew G. Knepley /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 130729f144ccSMatthew G. Knepley mpfr_set(psum, osum, MPFR_RNDN); 130829f144ccSMatthew G. Knepley mpfr_set(osum, sum, MPFR_RNDN); 130929f144ccSMatthew G. Knepley mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 131029f144ccSMatthew G. Knepley mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 131129f144ccSMatthew G. Knepley do { 131229f144ccSMatthew G. Knepley mpfr_set_si(kh, k, MPFR_RNDN); 131329f144ccSMatthew G. Knepley mpfr_mul(kh, kh, h, MPFR_RNDN); 131429f144ccSMatthew G. Knepley /* Weight */ 131529f144ccSMatthew G. Knepley mpfr_set(wk, h, MPFR_RNDN); 131629f144ccSMatthew G. Knepley mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 131729f144ccSMatthew G. Knepley mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 131829f144ccSMatthew G. Knepley mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 131929f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 132029f144ccSMatthew G. Knepley mpfr_sqr(tmp, tmp, MPFR_RNDN); 132129f144ccSMatthew G. Knepley mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 132229f144ccSMatthew G. Knepley mpfr_div(wk, wk, tmp, MPFR_RNDN); 132329f144ccSMatthew G. Knepley /* Abscissa */ 132429f144ccSMatthew G. Knepley mpfr_set_d(yk, 1.0, MPFR_RNDZ); 132529f144ccSMatthew G. Knepley mpfr_cosh(tmp, msinh, MPFR_RNDN); 132629f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 132729f144ccSMatthew G. Knepley mpfr_exp(tmp, msinh, MPFR_RNDN); 132829f144ccSMatthew G. Knepley mpfr_div(yk, yk, tmp, MPFR_RNDZ); 132929f144ccSMatthew G. Knepley /* Quadrature points */ 133029f144ccSMatthew G. Knepley mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 133129f144ccSMatthew G. Knepley mpfr_mul(lx, lx, alpha, MPFR_RNDU); 133229f144ccSMatthew G. Knepley mpfr_add(lx, lx, beta, MPFR_RNDU); 133329f144ccSMatthew G. Knepley mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 133429f144ccSMatthew G. Knepley mpfr_mul(rx, rx, alpha, MPFR_RNDD); 133529f144ccSMatthew G. Knepley mpfr_add(rx, rx, beta, MPFR_RNDD); 133629f144ccSMatthew G. Knepley /* Evaluation */ 133729f144ccSMatthew G. Knepley func(mpfr_get_d(lx, MPFR_RNDU), &lval); 133829f144ccSMatthew G. Knepley func(mpfr_get_d(rx, MPFR_RNDD), &rval); 133929f144ccSMatthew G. Knepley /* Update */ 134029f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 134129f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 134229f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 134329f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 134429f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 134529f144ccSMatthew G. Knepley mpfr_set(curTerm, tmp, MPFR_RNDN); 134629f144ccSMatthew G. Knepley mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 134729f144ccSMatthew G. Knepley mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 134829f144ccSMatthew G. Knepley mpfr_add(sum, sum, tmp, MPFR_RNDN); 134929f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 135029f144ccSMatthew G. Knepley mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 135129f144ccSMatthew G. Knepley mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 135229f144ccSMatthew G. Knepley ++k; 135329f144ccSMatthew G. Knepley /* Only need to evaluate every other point on refined levels */ 135429f144ccSMatthew G. Knepley if (l != 1) ++k; 135529f144ccSMatthew G. Knepley mpfr_log10(tmp, wk, MPFR_RNDN); 135629f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 1357c9f744b5SMatthew G. Knepley } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 135829f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, osum, MPFR_RNDN); 135929f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 136029f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 136129f144ccSMatthew G. Knepley d1 = mpfr_get_d(tmp, MPFR_RNDN); 136229f144ccSMatthew G. Knepley mpfr_sub(tmp, sum, psum, MPFR_RNDN); 136329f144ccSMatthew G. Knepley mpfr_abs(tmp, tmp, MPFR_RNDN); 136429f144ccSMatthew G. Knepley mpfr_log10(tmp, tmp, MPFR_RNDN); 136529f144ccSMatthew G. Knepley d2 = mpfr_get_d(tmp, MPFR_RNDN); 136629f144ccSMatthew G. Knepley mpfr_log10(tmp, maxTerm, MPFR_RNDN); 1367c9f744b5SMatthew G. Knepley d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 136829f144ccSMatthew G. Knepley mpfr_log10(tmp, curTerm, MPFR_RNDN); 136929f144ccSMatthew G. Knepley d4 = mpfr_get_d(tmp, MPFR_RNDN); 137029f144ccSMatthew G. Knepley d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1371b0649871SThomas Klotz } while (d < digits && l < 8); 137229f144ccSMatthew G. Knepley *sol = mpfr_get_d(sum, MPFR_RNDN); 137329f144ccSMatthew G. Knepley /* Cleanup */ 137429f144ccSMatthew G. Knepley mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 137529f144ccSMatthew G. Knepley PetscFunctionReturn(0); 137629f144ccSMatthew G. Knepley } 1377d525116cSMatthew G. Knepley #else 1378fbfcfee5SBarry Smith 1379d525116cSMatthew G. Knepley PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1380d525116cSMatthew G. Knepley { 1381d525116cSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 1382d525116cSMatthew G. Knepley } 138329f144ccSMatthew G. Knepley #endif 138429f144ccSMatthew G. Knepley 1385194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 1386194825f6SJed Brown * A in column-major format 1387194825f6SJed Brown * Ainv in row-major format 1388194825f6SJed Brown * tau has length m 1389194825f6SJed Brown * worksize must be >= max(1,n) 1390194825f6SJed Brown */ 1391194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1392194825f6SJed Brown { 1393194825f6SJed Brown PetscErrorCode ierr; 1394194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1395194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 1396194825f6SJed Brown 1397194825f6SJed Brown PetscFunctionBegin; 1398194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1399194825f6SJed Brown { 1400194825f6SJed Brown PetscInt i,j; 1401dcca6d9dSJed Brown ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 1402194825f6SJed Brown for (j=0; j<n; j++) { 1403194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 1404194825f6SJed Brown } 1405194825f6SJed Brown mstride = m; 1406194825f6SJed Brown } 1407194825f6SJed Brown #else 1408194825f6SJed Brown A = A_in; 1409194825f6SJed Brown Ainv = Ainv_out; 1410194825f6SJed Brown #endif 1411194825f6SJed Brown 1412194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 1413194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 1414194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 1415194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 1416194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1417001a771dSBarry Smith PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1418194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 1419194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1420194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1421194825f6SJed Brown 1422194825f6SJed Brown /* Extract an explicit representation of Q */ 1423194825f6SJed Brown Q = Ainv; 1424580bdb30SBarry Smith ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr); 1425194825f6SJed Brown K = N; /* full rank */ 1426c964aadfSJose E. Roman PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1427194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1428194825f6SJed Brown 1429194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1430194825f6SJed Brown Alpha = 1.0; 1431194825f6SJed Brown ldb = lda; 1432001a771dSBarry Smith PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 1433194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 1434194825f6SJed Brown 1435194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 1436194825f6SJed Brown { 1437194825f6SJed Brown PetscInt i; 1438194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 1439194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 1440194825f6SJed Brown } 1441194825f6SJed Brown #endif 1442194825f6SJed Brown PetscFunctionReturn(0); 1443194825f6SJed Brown } 1444194825f6SJed Brown 1445194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 1446194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 1447194825f6SJed Brown { 1448194825f6SJed Brown PetscErrorCode ierr; 1449194825f6SJed Brown PetscReal *Bv; 1450194825f6SJed Brown PetscInt i,j; 1451194825f6SJed Brown 1452194825f6SJed Brown PetscFunctionBegin; 1453785e854fSJed Brown ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 1454194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 1455194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 1456194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 1457194825f6SJed Brown for (i=0; i<ninterval; i++) { 1458194825f6SJed Brown for (j=0; j<ndegree; j++) { 1459194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1460194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1461194825f6SJed Brown } 1462194825f6SJed Brown } 1463194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 1464194825f6SJed Brown PetscFunctionReturn(0); 1465194825f6SJed Brown } 1466194825f6SJed Brown 1467194825f6SJed Brown /*@ 1468194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 1469194825f6SJed Brown 1470194825f6SJed Brown Not Collective 1471194825f6SJed Brown 1472194825f6SJed Brown Input Arguments: 1473194825f6SJed Brown + degree - degree of reconstruction polynomial 1474194825f6SJed Brown . nsource - number of source intervals 1475194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 1476194825f6SJed Brown . ntarget - number of target intervals 1477194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 1478194825f6SJed Brown 1479194825f6SJed Brown Output Arguments: 1480194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 1481194825f6SJed Brown 1482194825f6SJed Brown Level: advanced 1483194825f6SJed Brown 1484194825f6SJed Brown .seealso: PetscDTLegendreEval() 1485194825f6SJed Brown @*/ 1486194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 1487194825f6SJed Brown { 1488194825f6SJed Brown PetscErrorCode ierr; 1489194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 1490194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 1491194825f6SJed Brown PetscScalar *tau,*work; 1492194825f6SJed Brown 1493194825f6SJed Brown PetscFunctionBegin; 1494194825f6SJed Brown PetscValidRealPointer(sourcex,3); 1495194825f6SJed Brown PetscValidRealPointer(targetx,5); 1496194825f6SJed Brown PetscValidRealPointer(R,6); 1497194825f6SJed 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); 1498194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 1499194825f6SJed Brown for (i=0; i<nsource; i++) { 150057622a8eSBarry 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]); 1501194825f6SJed Brown } 1502194825f6SJed Brown for (i=0; i<ntarget; i++) { 150357622a8eSBarry 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]); 1504194825f6SJed Brown } 1505194825f6SJed Brown #endif 1506194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 1507194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 1508194825f6SJed Brown center = (xmin + xmax)/2; 1509194825f6SJed Brown hscale = (xmax - xmin)/2; 1510194825f6SJed Brown worksize = nsource; 1511dcca6d9dSJed Brown ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 1512dcca6d9dSJed Brown ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 1513194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 1514194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 1515194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 1516194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 1517194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 1518194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 1519194825f6SJed Brown for (i=0; i<ntarget; i++) { 1520194825f6SJed Brown PetscReal rowsum = 0; 1521194825f6SJed Brown for (j=0; j<nsource; j++) { 1522194825f6SJed Brown PetscReal sum = 0; 1523194825f6SJed Brown for (k=0; k<degree+1; k++) { 1524194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 1525194825f6SJed Brown } 1526194825f6SJed Brown R[i*nsource+j] = sum; 1527194825f6SJed Brown rowsum += sum; 1528194825f6SJed Brown } 1529194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 1530194825f6SJed Brown } 1531194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 1532194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 1533194825f6SJed Brown PetscFunctionReturn(0); 1534194825f6SJed Brown } 1535916e780bShannah_mairs 1536916e780bShannah_mairs /*@C 1537916e780bShannah_mairs PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points 1538916e780bShannah_mairs 1539916e780bShannah_mairs Not Collective 1540916e780bShannah_mairs 1541916e780bShannah_mairs Input Parameter: 1542916e780bShannah_mairs + n - the number of GLL nodes 1543916e780bShannah_mairs . nodes - the GLL nodes 1544916e780bShannah_mairs . weights - the GLL weights 1545916e780bShannah_mairs . f - the function values at the nodes 1546916e780bShannah_mairs 1547916e780bShannah_mairs Output Parameter: 1548916e780bShannah_mairs . in - the value of the integral 1549916e780bShannah_mairs 1550916e780bShannah_mairs Level: beginner 1551916e780bShannah_mairs 1552916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature() 1553916e780bShannah_mairs 1554916e780bShannah_mairs @*/ 1555916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in) 1556916e780bShannah_mairs { 1557916e780bShannah_mairs PetscInt i; 1558916e780bShannah_mairs 1559916e780bShannah_mairs PetscFunctionBegin; 1560916e780bShannah_mairs *in = 0.; 1561916e780bShannah_mairs for (i=0; i<n; i++) { 1562916e780bShannah_mairs *in += f[i]*f[i]*weights[i]; 1563916e780bShannah_mairs } 1564916e780bShannah_mairs PetscFunctionReturn(0); 1565916e780bShannah_mairs } 1566916e780bShannah_mairs 1567916e780bShannah_mairs /*@C 1568916e780bShannah_mairs PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element 1569916e780bShannah_mairs 1570916e780bShannah_mairs Not Collective 1571916e780bShannah_mairs 1572916e780bShannah_mairs Input Parameter: 1573916e780bShannah_mairs + n - the number of GLL nodes 1574916e780bShannah_mairs . nodes - the GLL nodes 1575916e780bShannah_mairs . weights - the GLL weights 1576916e780bShannah_mairs 1577916e780bShannah_mairs Output Parameter: 1578916e780bShannah_mairs . A - the stiffness element 1579916e780bShannah_mairs 1580916e780bShannah_mairs Level: beginner 1581916e780bShannah_mairs 1582916e780bShannah_mairs Notes: 1583916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy() 1584916e780bShannah_mairs 1585916e780bShannah_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) 1586916e780bShannah_mairs 1587916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1588916e780bShannah_mairs 1589916e780bShannah_mairs @*/ 1590916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1591916e780bShannah_mairs { 1592916e780bShannah_mairs PetscReal **A; 1593916e780bShannah_mairs PetscErrorCode ierr; 1594916e780bShannah_mairs const PetscReal *gllnodes = nodes; 1595916e780bShannah_mairs const PetscInt p = n-1; 1596916e780bShannah_mairs PetscReal z0,z1,z2 = -1,x,Lpj,Lpr; 1597916e780bShannah_mairs PetscInt i,j,nn,r; 1598916e780bShannah_mairs 1599916e780bShannah_mairs PetscFunctionBegin; 1600916e780bShannah_mairs ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1601916e780bShannah_mairs ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1602916e780bShannah_mairs for (i=1; i<n; i++) A[i] = A[i-1]+n; 1603916e780bShannah_mairs 1604916e780bShannah_mairs for (j=1; j<p; j++) { 1605916e780bShannah_mairs x = gllnodes[j]; 1606916e780bShannah_mairs z0 = 1.; 1607916e780bShannah_mairs z1 = x; 1608916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1609916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1610916e780bShannah_mairs z0 = z1; 1611916e780bShannah_mairs z1 = z2; 1612916e780bShannah_mairs } 1613916e780bShannah_mairs Lpj=z2; 1614916e780bShannah_mairs for (r=1; r<p; r++) { 1615916e780bShannah_mairs if (r == j) { 1616916e780bShannah_mairs A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj); 1617916e780bShannah_mairs } else { 1618916e780bShannah_mairs x = gllnodes[r]; 1619916e780bShannah_mairs z0 = 1.; 1620916e780bShannah_mairs z1 = x; 1621916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1622916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1623916e780bShannah_mairs z0 = z1; 1624916e780bShannah_mairs z1 = z2; 1625916e780bShannah_mairs } 1626916e780bShannah_mairs Lpr = z2; 1627916e780bShannah_mairs A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r])); 1628916e780bShannah_mairs } 1629916e780bShannah_mairs } 1630916e780bShannah_mairs } 1631916e780bShannah_mairs for (j=1; j<p+1; j++) { 1632916e780bShannah_mairs x = gllnodes[j]; 1633916e780bShannah_mairs z0 = 1.; 1634916e780bShannah_mairs z1 = x; 1635916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1636916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1637916e780bShannah_mairs z0 = z1; 1638916e780bShannah_mairs z1 = z2; 1639916e780bShannah_mairs } 1640916e780bShannah_mairs Lpj = z2; 1641916e780bShannah_mairs A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j])); 1642916e780bShannah_mairs A[0][j] = A[j][0]; 1643916e780bShannah_mairs } 1644916e780bShannah_mairs for (j=0; j<p; j++) { 1645916e780bShannah_mairs x = gllnodes[j]; 1646916e780bShannah_mairs z0 = 1.; 1647916e780bShannah_mairs z1 = x; 1648916e780bShannah_mairs for (nn=1; nn<p; nn++) { 1649916e780bShannah_mairs z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1650916e780bShannah_mairs z0 = z1; 1651916e780bShannah_mairs z1 = z2; 1652916e780bShannah_mairs } 1653916e780bShannah_mairs Lpj=z2; 1654916e780bShannah_mairs 1655916e780bShannah_mairs A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j])); 1656916e780bShannah_mairs A[j][p] = A[p][j]; 1657916e780bShannah_mairs } 1658916e780bShannah_mairs A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.; 1659916e780bShannah_mairs A[p][p]=A[0][0]; 1660916e780bShannah_mairs *AA = A; 1661916e780bShannah_mairs PetscFunctionReturn(0); 1662916e780bShannah_mairs } 1663916e780bShannah_mairs 1664916e780bShannah_mairs /*@C 1665916e780bShannah_mairs PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element 1666916e780bShannah_mairs 1667916e780bShannah_mairs Not Collective 1668916e780bShannah_mairs 1669916e780bShannah_mairs Input Parameter: 1670916e780bShannah_mairs + n - the number of GLL nodes 1671916e780bShannah_mairs . nodes - the GLL nodes 1672916e780bShannah_mairs . weights - the GLL weightss 1673916e780bShannah_mairs - A - the stiffness element 1674916e780bShannah_mairs 1675916e780bShannah_mairs Level: beginner 1676916e780bShannah_mairs 1677916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate() 1678916e780bShannah_mairs 1679916e780bShannah_mairs @*/ 1680916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1681916e780bShannah_mairs { 1682916e780bShannah_mairs PetscErrorCode ierr; 1683916e780bShannah_mairs 1684916e780bShannah_mairs PetscFunctionBegin; 1685916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1686916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1687916e780bShannah_mairs *AA = NULL; 1688916e780bShannah_mairs PetscFunctionReturn(0); 1689916e780bShannah_mairs } 1690916e780bShannah_mairs 1691916e780bShannah_mairs /*@C 1692916e780bShannah_mairs PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element 1693916e780bShannah_mairs 1694916e780bShannah_mairs Not Collective 1695916e780bShannah_mairs 1696916e780bShannah_mairs Input Parameter: 1697916e780bShannah_mairs + n - the number of GLL nodes 1698916e780bShannah_mairs . nodes - the GLL nodes 1699916e780bShannah_mairs . weights - the GLL weights 1700916e780bShannah_mairs 1701916e780bShannah_mairs Output Parameter: 1702916e780bShannah_mairs . AA - the stiffness element 1703916e780bShannah_mairs - AAT - the transpose of AA (pass in NULL if you do not need this array) 1704916e780bShannah_mairs 1705916e780bShannah_mairs Level: beginner 1706916e780bShannah_mairs 1707916e780bShannah_mairs Notes: 1708916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementGradientDestroy() 1709916e780bShannah_mairs 1710916e780bShannah_mairs You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1711916e780bShannah_mairs 1712916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1713916e780bShannah_mairs 1714916e780bShannah_mairs @*/ 1715916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1716916e780bShannah_mairs { 1717916e780bShannah_mairs PetscReal **A, **AT = NULL; 1718916e780bShannah_mairs PetscErrorCode ierr; 1719916e780bShannah_mairs const PetscReal *gllnodes = nodes; 1720916e780bShannah_mairs const PetscInt p = n-1; 1721916e780bShannah_mairs PetscReal q,qp,Li, Lj,d0; 1722916e780bShannah_mairs PetscInt i,j; 1723916e780bShannah_mairs 1724916e780bShannah_mairs PetscFunctionBegin; 1725916e780bShannah_mairs ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1726916e780bShannah_mairs ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1727916e780bShannah_mairs for (i=1; i<n; i++) A[i] = A[i-1]+n; 1728916e780bShannah_mairs 1729916e780bShannah_mairs if (AAT) { 1730916e780bShannah_mairs ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr); 1731916e780bShannah_mairs ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr); 1732916e780bShannah_mairs for (i=1; i<n; i++) AT[i] = AT[i-1]+n; 1733916e780bShannah_mairs } 1734916e780bShannah_mairs 1735916e780bShannah_mairs if (n==1) {A[0][0] = 0.;} 1736916e780bShannah_mairs d0 = (PetscReal)p*((PetscReal)p+1.)/4.; 1737916e780bShannah_mairs for (i=0; i<n; i++) { 1738916e780bShannah_mairs for (j=0; j<n; j++) { 1739916e780bShannah_mairs A[i][j] = 0.; 1740fdd31e58Shannah_mairs qAndLEvaluation(p,gllnodes[i],&q,&qp,&Li); 1741fdd31e58Shannah_mairs qAndLEvaluation(p,gllnodes[j],&q,&qp,&Lj); 1742916e780bShannah_mairs if (i!=j) A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j])); 1743916e780bShannah_mairs if ((j==i) && (i==0)) A[i][j] = -d0; 1744916e780bShannah_mairs if (j==i && i==p) A[i][j] = d0; 1745916e780bShannah_mairs if (AT) AT[j][i] = A[i][j]; 1746916e780bShannah_mairs } 1747916e780bShannah_mairs } 1748916e780bShannah_mairs if (AAT) *AAT = AT; 1749916e780bShannah_mairs *AA = A; 1750916e780bShannah_mairs PetscFunctionReturn(0); 1751916e780bShannah_mairs } 1752916e780bShannah_mairs 1753916e780bShannah_mairs /*@C 1754916e780bShannah_mairs PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate() 1755916e780bShannah_mairs 1756916e780bShannah_mairs Not Collective 1757916e780bShannah_mairs 1758916e780bShannah_mairs Input Parameter: 1759916e780bShannah_mairs + n - the number of GLL nodes 1760916e780bShannah_mairs . nodes - the GLL nodes 1761916e780bShannah_mairs . weights - the GLL weights 1762916e780bShannah_mairs . AA - the stiffness element 1763916e780bShannah_mairs - AAT - the transpose of the element 1764916e780bShannah_mairs 1765916e780bShannah_mairs Level: beginner 1766916e780bShannah_mairs 1767916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate() 1768916e780bShannah_mairs 1769916e780bShannah_mairs @*/ 1770916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1771916e780bShannah_mairs { 1772916e780bShannah_mairs PetscErrorCode ierr; 1773916e780bShannah_mairs 1774916e780bShannah_mairs PetscFunctionBegin; 1775916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1776916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1777916e780bShannah_mairs *AA = NULL; 1778916e780bShannah_mairs if (*AAT) { 1779916e780bShannah_mairs ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr); 1780916e780bShannah_mairs ierr = PetscFree(*AAT);CHKERRQ(ierr); 1781916e780bShannah_mairs *AAT = NULL; 1782916e780bShannah_mairs } 1783916e780bShannah_mairs PetscFunctionReturn(0); 1784916e780bShannah_mairs } 1785916e780bShannah_mairs 1786916e780bShannah_mairs /*@C 1787916e780bShannah_mairs PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element 1788916e780bShannah_mairs 1789916e780bShannah_mairs Not Collective 1790916e780bShannah_mairs 1791916e780bShannah_mairs Input Parameter: 1792916e780bShannah_mairs + n - the number of GLL nodes 1793916e780bShannah_mairs . nodes - the GLL nodes 1794916e780bShannah_mairs . weights - the GLL weightss 1795916e780bShannah_mairs 1796916e780bShannah_mairs Output Parameter: 1797916e780bShannah_mairs . AA - the stiffness element 1798916e780bShannah_mairs 1799916e780bShannah_mairs Level: beginner 1800916e780bShannah_mairs 1801916e780bShannah_mairs Notes: 1802916e780bShannah_mairs Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy() 1803916e780bShannah_mairs 1804916e780bShannah_mairs This is the same as the Gradient operator multiplied by the diagonal mass matrix 1805916e780bShannah_mairs 1806916e780bShannah_mairs You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1807916e780bShannah_mairs 1808916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy() 1809916e780bShannah_mairs 1810916e780bShannah_mairs @*/ 1811916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1812916e780bShannah_mairs { 1813916e780bShannah_mairs PetscReal **D; 1814916e780bShannah_mairs PetscErrorCode ierr; 1815916e780bShannah_mairs const PetscReal *gllweights = weights; 1816916e780bShannah_mairs const PetscInt glln = n; 1817916e780bShannah_mairs PetscInt i,j; 1818916e780bShannah_mairs 1819916e780bShannah_mairs PetscFunctionBegin; 1820916e780bShannah_mairs ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr); 1821916e780bShannah_mairs for (i=0; i<glln; i++){ 1822916e780bShannah_mairs for (j=0; j<glln; j++) { 1823916e780bShannah_mairs D[i][j] = gllweights[i]*D[i][j]; 1824916e780bShannah_mairs } 1825916e780bShannah_mairs } 1826916e780bShannah_mairs *AA = D; 1827916e780bShannah_mairs PetscFunctionReturn(0); 1828916e780bShannah_mairs } 1829916e780bShannah_mairs 1830916e780bShannah_mairs /*@C 1831916e780bShannah_mairs PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element 1832916e780bShannah_mairs 1833916e780bShannah_mairs Not Collective 1834916e780bShannah_mairs 1835916e780bShannah_mairs Input Parameter: 1836916e780bShannah_mairs + n - the number of GLL nodes 1837916e780bShannah_mairs . nodes - the GLL nodes 1838916e780bShannah_mairs . weights - the GLL weights 1839916e780bShannah_mairs - A - advection 1840916e780bShannah_mairs 1841916e780bShannah_mairs Level: beginner 1842916e780bShannah_mairs 1843916e780bShannah_mairs .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate() 1844916e780bShannah_mairs 1845916e780bShannah_mairs @*/ 1846916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1847916e780bShannah_mairs { 1848916e780bShannah_mairs PetscErrorCode ierr; 1849916e780bShannah_mairs 1850916e780bShannah_mairs PetscFunctionBegin; 1851916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1852916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1853916e780bShannah_mairs *AA = NULL; 1854916e780bShannah_mairs PetscFunctionReturn(0); 1855916e780bShannah_mairs } 1856916e780bShannah_mairs 1857916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1858916e780bShannah_mairs { 1859916e780bShannah_mairs PetscReal **A; 1860916e780bShannah_mairs PetscErrorCode ierr; 1861916e780bShannah_mairs const PetscReal *gllweights = weights; 1862916e780bShannah_mairs const PetscInt glln = n; 1863916e780bShannah_mairs PetscInt i,j; 1864916e780bShannah_mairs 1865916e780bShannah_mairs PetscFunctionBegin; 1866916e780bShannah_mairs ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr); 1867916e780bShannah_mairs ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr); 1868916e780bShannah_mairs for (i=1; i<glln; i++) A[i] = A[i-1]+glln; 1869916e780bShannah_mairs if (glln==1) {A[0][0] = 0.;} 1870916e780bShannah_mairs for (i=0; i<glln; i++) { 1871916e780bShannah_mairs for (j=0; j<glln; j++) { 1872916e780bShannah_mairs A[i][j] = 0.; 1873916e780bShannah_mairs if (j==i) A[i][j] = gllweights[i]; 1874916e780bShannah_mairs } 1875916e780bShannah_mairs } 1876916e780bShannah_mairs *AA = A; 1877916e780bShannah_mairs PetscFunctionReturn(0); 1878916e780bShannah_mairs } 1879916e780bShannah_mairs 1880916e780bShannah_mairs PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1881916e780bShannah_mairs { 1882916e780bShannah_mairs PetscErrorCode ierr; 1883916e780bShannah_mairs 1884916e780bShannah_mairs PetscFunctionBegin; 1885916e780bShannah_mairs ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1886916e780bShannah_mairs ierr = PetscFree(*AA);CHKERRQ(ierr); 1887916e780bShannah_mairs *AA = NULL; 1888916e780bShannah_mairs PetscFunctionReturn(0); 1889916e780bShannah_mairs } 1890916e780bShannah_mairs 1891