1 /* Discretization tools */ 2 3 #include <petscdt.h> /*I "petscdt.h" I*/ 4 #include <petscblaslapack.h> 5 #include <petsc/private/petscimpl.h> 6 #include <petsc/private/dtimpl.h> 7 #include <petscviewer.h> 8 #include <petscdmplex.h> 9 #include <petscdmshell.h> 10 11 #if defined(PETSC_HAVE_MPFR) 12 #include <mpfr.h> 13 #endif 14 15 static PetscBool GolubWelschCite = PETSC_FALSE; 16 const char GolubWelschCitation[] = "@article{GolubWelsch1969,\n" 17 " author = {Golub and Welsch},\n" 18 " title = {Calculation of Quadrature Rules},\n" 19 " journal = {Math. Comp.},\n" 20 " volume = {23},\n" 21 " number = {106},\n" 22 " pages = {221--230},\n" 23 " year = {1969}\n}\n"; 24 25 static PetscBool gaussQuadratureNewton = PETSC_FALSE; 26 27 28 PetscClassId PETSCQUADRATURE_CLASSID = 0; 29 30 /*@ 31 PetscQuadratureCreate - Create a PetscQuadrature object 32 33 Collective 34 35 Input Parameter: 36 . comm - The communicator for the PetscQuadrature object 37 38 Output Parameter: 39 . q - The PetscQuadrature object 40 41 Level: beginner 42 43 .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 44 @*/ 45 PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 46 { 47 PetscErrorCode ierr; 48 49 PetscFunctionBegin; 50 PetscValidPointer(q, 2); 51 ierr = DMInitializePackage();CHKERRQ(ierr); 52 ierr = PetscHeaderCreate(*q,PETSCQUADRATURE_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 53 (*q)->dim = -1; 54 (*q)->Nc = 1; 55 (*q)->order = -1; 56 (*q)->numPoints = 0; 57 (*q)->points = NULL; 58 (*q)->weights = NULL; 59 PetscFunctionReturn(0); 60 } 61 62 /*@ 63 PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 64 65 Collective on q 66 67 Input Parameter: 68 . q - The PetscQuadrature object 69 70 Output Parameter: 71 . r - The new PetscQuadrature object 72 73 Level: beginner 74 75 .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 76 @*/ 77 PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 78 { 79 PetscInt order, dim, Nc, Nq; 80 const PetscReal *points, *weights; 81 PetscReal *p, *w; 82 PetscErrorCode ierr; 83 84 PetscFunctionBegin; 85 PetscValidPointer(q, 2); 86 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 87 ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 88 ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 89 ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 90 ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 91 ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr); 92 ierr = PetscArraycpy(p, points, Nq*dim);CHKERRQ(ierr); 93 ierr = PetscArraycpy(w, weights, Nc * Nq);CHKERRQ(ierr); 94 ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr); 95 PetscFunctionReturn(0); 96 } 97 98 /*@ 99 PetscQuadratureDestroy - Destroys a PetscQuadrature object 100 101 Collective on q 102 103 Input Parameter: 104 . q - The PetscQuadrature object 105 106 Level: beginner 107 108 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 109 @*/ 110 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 111 { 112 PetscErrorCode ierr; 113 114 PetscFunctionBegin; 115 if (!*q) PetscFunctionReturn(0); 116 PetscValidHeaderSpecific((*q),PETSCQUADRATURE_CLASSID,1); 117 if (--((PetscObject)(*q))->refct > 0) { 118 *q = NULL; 119 PetscFunctionReturn(0); 120 } 121 ierr = PetscFree((*q)->points);CHKERRQ(ierr); 122 ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 123 ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 124 PetscFunctionReturn(0); 125 } 126 127 /*@ 128 PetscQuadratureGetOrder - Return the order of the method 129 130 Not collective 131 132 Input Parameter: 133 . q - The PetscQuadrature object 134 135 Output Parameter: 136 . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 137 138 Level: intermediate 139 140 .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 141 @*/ 142 PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 143 { 144 PetscFunctionBegin; 145 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 146 PetscValidPointer(order, 2); 147 *order = q->order; 148 PetscFunctionReturn(0); 149 } 150 151 /*@ 152 PetscQuadratureSetOrder - Return the order of the method 153 154 Not collective 155 156 Input Parameters: 157 + q - The PetscQuadrature object 158 - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 159 160 Level: intermediate 161 162 .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 163 @*/ 164 PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 165 { 166 PetscFunctionBegin; 167 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 168 q->order = order; 169 PetscFunctionReturn(0); 170 } 171 172 /*@ 173 PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated 174 175 Not collective 176 177 Input Parameter: 178 . q - The PetscQuadrature object 179 180 Output Parameter: 181 . Nc - The number of components 182 183 Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 184 185 Level: intermediate 186 187 .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 188 @*/ 189 PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc) 190 { 191 PetscFunctionBegin; 192 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 193 PetscValidPointer(Nc, 2); 194 *Nc = q->Nc; 195 PetscFunctionReturn(0); 196 } 197 198 /*@ 199 PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated 200 201 Not collective 202 203 Input Parameters: 204 + q - The PetscQuadrature object 205 - Nc - The number of components 206 207 Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 208 209 Level: intermediate 210 211 .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 212 @*/ 213 PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc) 214 { 215 PetscFunctionBegin; 216 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 217 q->Nc = Nc; 218 PetscFunctionReturn(0); 219 } 220 221 /*@C 222 PetscQuadratureGetData - Returns the data defining the quadrature 223 224 Not collective 225 226 Input Parameter: 227 . q - The PetscQuadrature object 228 229 Output Parameters: 230 + dim - The spatial dimension 231 . Nc - The number of components 232 . npoints - The number of quadrature points 233 . points - The coordinates of each quadrature point 234 - weights - The weight of each quadrature point 235 236 Level: intermediate 237 238 Fortran Notes: 239 From Fortran you must call PetscQuadratureRestoreData() when you are done with the data 240 241 .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 242 @*/ 243 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 244 { 245 PetscFunctionBegin; 246 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 247 if (dim) { 248 PetscValidPointer(dim, 2); 249 *dim = q->dim; 250 } 251 if (Nc) { 252 PetscValidPointer(Nc, 3); 253 *Nc = q->Nc; 254 } 255 if (npoints) { 256 PetscValidPointer(npoints, 4); 257 *npoints = q->numPoints; 258 } 259 if (points) { 260 PetscValidPointer(points, 5); 261 *points = q->points; 262 } 263 if (weights) { 264 PetscValidPointer(weights, 6); 265 *weights = q->weights; 266 } 267 PetscFunctionReturn(0); 268 } 269 270 static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[]) 271 { 272 PetscScalar *Js, *Jinvs; 273 PetscInt i, j, k; 274 PetscBLASInt bm, bn, info; 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 ierr = PetscBLASIntCast(m, &bm);CHKERRQ(ierr); 279 ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr); 280 #if defined(PETSC_USE_COMPLEX) 281 ierr = PetscMalloc2(m*n, &Js, m*n, &Jinvs);CHKERRQ(ierr); 282 for (i = 0; i < m*n; i++) Js[i] = J[i]; 283 #else 284 Js = (PetscReal *) J; 285 Jinvs = Jinv; 286 #endif 287 if (m == n) { 288 PetscBLASInt *pivots; 289 PetscScalar *W; 290 291 ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr); 292 293 ierr = PetscArraycpy(Jinvs, Js, m * m);CHKERRQ(ierr); 294 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info)); 295 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 296 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info)); 297 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 298 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 299 } else if (m < n) { 300 PetscScalar *JJT; 301 PetscBLASInt *pivots; 302 PetscScalar *W; 303 304 ierr = PetscMalloc1(m*m, &JJT);CHKERRQ(ierr); 305 ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr); 306 for (i = 0; i < m; i++) { 307 for (j = 0; j < m; j++) { 308 PetscScalar val = 0.; 309 310 for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k]; 311 JJT[i * m + j] = val; 312 } 313 } 314 315 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info)); 316 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 317 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info)); 318 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 319 for (i = 0; i < n; i++) { 320 for (j = 0; j < m; j++) { 321 PetscScalar val = 0.; 322 323 for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j]; 324 Jinvs[i * m + j] = val; 325 } 326 } 327 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 328 ierr = PetscFree(JJT);CHKERRQ(ierr); 329 } else { 330 PetscScalar *JTJ; 331 PetscBLASInt *pivots; 332 PetscScalar *W; 333 334 ierr = PetscMalloc1(n*n, &JTJ);CHKERRQ(ierr); 335 ierr = PetscMalloc2(n, &pivots, n, &W);CHKERRQ(ierr); 336 for (i = 0; i < n; i++) { 337 for (j = 0; j < n; j++) { 338 PetscScalar val = 0.; 339 340 for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j]; 341 JTJ[i * n + j] = val; 342 } 343 } 344 345 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bm, pivots, &info)); 346 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 347 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info)); 348 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 349 for (i = 0; i < n; i++) { 350 for (j = 0; j < m; j++) { 351 PetscScalar val = 0.; 352 353 for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k]; 354 Jinvs[i * m + j] = val; 355 } 356 } 357 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 358 ierr = PetscFree(JTJ);CHKERRQ(ierr); 359 } 360 #if defined(PETSC_USE_COMPLEX) 361 for (i = 0; i < m*n; i++) Jinv[i] = PetscRealPart(Jinvs[i]); 362 ierr = PetscFree2(Js, Jinvs);CHKERRQ(ierr); 363 #endif 364 PetscFunctionReturn(0); 365 } 366 367 /*@ 368 PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation. 369 370 Collecive on PetscQuadrature 371 372 Input Arguments: 373 + q - the quadrature functional 374 . imageDim - the dimension of the image of the transformation 375 . origin - a point in the original space 376 . originImage - the image of the origin under the transformation 377 . J - the Jacobian of the image: an [imageDim x dim] matrix in row major order 378 - formDegree - transform the quadrature weights as k-forms of this form degree (if the number of components is a multiple of (dim choose formDegree), it is assumed that they represent multiple k-forms) [see PetscDTAltVPullback() for interpretation of formDegree] 379 380 Output Arguments: 381 . 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. 382 383 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. 384 385 .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 386 @*/ 387 PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formDegree, PetscQuadrature *Jinvstarq) 388 { 389 PetscInt dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c; 390 const PetscReal *points; 391 const PetscReal *weights; 392 PetscReal *imagePoints, *imageWeights; 393 PetscReal *Jinv; 394 PetscReal *Jinvstar; 395 PetscErrorCode ierr; 396 397 PetscFunctionBegin; 398 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 399 if (imageDim < PetscAbsInt(formDegree)) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %D-form in %D dimensions", PetscAbsInt(formDegree), imageDim); 400 ierr = PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights);CHKERRQ(ierr); 401 ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &formSize);CHKERRQ(ierr); 402 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); 403 Ncopies = Nc / formSize; 404 ierr = PetscDTBinomialInt(imageDim, PetscAbsInt(formDegree), &imageFormSize);CHKERRQ(ierr); 405 imageNc = Ncopies * imageFormSize; 406 ierr = PetscMalloc1(Npoints * imageDim, &imagePoints);CHKERRQ(ierr); 407 ierr = PetscMalloc1(Npoints * imageNc, &imageWeights);CHKERRQ(ierr); 408 ierr = PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar);CHKERRQ(ierr); 409 ierr = PetscDTJacobianInverse_Internal(dim, imageDim, J, Jinv);CHKERRQ(ierr); 410 ierr = PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formDegree, Jinvstar);CHKERRQ(ierr); 411 for (pt = 0; pt < Npoints; pt++) { 412 const PetscReal *point = &points[pt * dim]; 413 PetscReal *imagePoint = &imagePoints[pt * imageDim]; 414 415 for (i = 0; i < imageDim; i++) { 416 PetscReal val = originImage[i]; 417 418 for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]); 419 imagePoint[i] = val; 420 } 421 for (c = 0; c < Ncopies; c++) { 422 const PetscReal *form = &weights[pt * Nc + c * formSize]; 423 PetscReal *imageForm = &imageWeights[pt * imageNc + c * imageFormSize]; 424 425 for (i = 0; i < imageFormSize; i++) { 426 PetscReal val = 0.; 427 428 for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j]; 429 imageForm[i] = val; 430 } 431 } 432 } 433 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq);CHKERRQ(ierr); 434 ierr = PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights);CHKERRQ(ierr); 435 ierr = PetscFree2(Jinv, Jinvstar);CHKERRQ(ierr); 436 PetscFunctionReturn(0); 437 } 438 439 /*@C 440 PetscQuadratureSetData - Sets the data defining the quadrature 441 442 Not collective 443 444 Input Parameters: 445 + q - The PetscQuadrature object 446 . dim - The spatial dimension 447 . Nc - The number of components 448 . npoints - The number of quadrature points 449 . points - The coordinates of each quadrature point 450 - weights - The weight of each quadrature point 451 452 Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them. 453 454 Level: intermediate 455 456 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 457 @*/ 458 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 459 { 460 PetscFunctionBegin; 461 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 462 if (dim >= 0) q->dim = dim; 463 if (Nc >= 0) q->Nc = Nc; 464 if (npoints >= 0) q->numPoints = npoints; 465 if (points) { 466 PetscValidPointer(points, 4); 467 q->points = points; 468 } 469 if (weights) { 470 PetscValidPointer(weights, 5); 471 q->weights = weights; 472 } 473 PetscFunctionReturn(0); 474 } 475 476 static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v) 477 { 478 PetscInt q, d, c; 479 PetscViewerFormat format; 480 PetscErrorCode ierr; 481 482 PetscFunctionBegin; 483 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);} 484 else {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);} 485 ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 486 if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0); 487 for (q = 0; q < quad->numPoints; ++q) { 488 ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr); 489 ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr); 490 for (d = 0; d < quad->dim; ++d) { 491 if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 492 ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 493 } 494 ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr); 495 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);} 496 for (c = 0; c < quad->Nc; ++c) { 497 if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 498 ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr); 499 } 500 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);} 501 ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr); 502 ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr); 503 } 504 PetscFunctionReturn(0); 505 } 506 507 /*@C 508 PetscQuadratureView - Views a PetscQuadrature object 509 510 Collective on quad 511 512 Input Parameters: 513 + quad - The PetscQuadrature object 514 - viewer - The PetscViewer object 515 516 Level: beginner 517 518 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 519 @*/ 520 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 521 { 522 PetscBool iascii; 523 PetscErrorCode ierr; 524 525 PetscFunctionBegin; 526 PetscValidHeader(quad, 1); 527 if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 528 if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);} 529 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 530 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 531 if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);} 532 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 533 PetscFunctionReturn(0); 534 } 535 536 /*@C 537 PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 538 539 Not collective 540 541 Input Parameter: 542 + q - The original PetscQuadrature 543 . numSubelements - The number of subelements the original element is divided into 544 . v0 - An array of the initial points for each subelement 545 - jac - An array of the Jacobian mappings from the reference to each subelement 546 547 Output Parameters: 548 . dim - The dimension 549 550 Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 551 552 Not available from Fortran 553 554 Level: intermediate 555 556 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 557 @*/ 558 PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 559 { 560 const PetscReal *points, *weights; 561 PetscReal *pointsRef, *weightsRef; 562 PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e; 563 PetscErrorCode ierr; 564 565 PetscFunctionBegin; 566 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 567 PetscValidPointer(v0, 3); 568 PetscValidPointer(jac, 4); 569 PetscValidPointer(qref, 5); 570 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 571 ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 572 ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 573 npointsRef = npoints*numSubelements; 574 ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 575 ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr); 576 for (c = 0; c < numSubelements; ++c) { 577 for (p = 0; p < npoints; ++p) { 578 for (d = 0; d < dim; ++d) { 579 pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 580 for (e = 0; e < dim; ++e) { 581 pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 582 } 583 } 584 /* Could also use detJ here */ 585 for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements; 586 } 587 } 588 ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 589 ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 590 PetscFunctionReturn(0); 591 } 592 593 /*@ 594 PetscDTLegendreEval - evaluate Legendre polynomial at points 595 596 Not Collective 597 598 Input Arguments: 599 + npoints - number of spatial points to evaluate at 600 . points - array of locations to evaluate at 601 . ndegree - number of basis degrees to evaluate 602 - degrees - sorted array of degrees to evaluate 603 604 Output Arguments: 605 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 606 . D - row-oriented derivative evaluation matrix (or NULL) 607 - D2 - row-oriented second derivative evaluation matrix (or NULL) 608 609 Level: intermediate 610 611 .seealso: PetscDTGaussQuadrature() 612 @*/ 613 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 614 { 615 PetscInt i,maxdegree; 616 617 PetscFunctionBegin; 618 if (!npoints || !ndegree) PetscFunctionReturn(0); 619 maxdegree = degrees[ndegree-1]; 620 for (i=0; i<npoints; i++) { 621 PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 622 PetscInt j,k; 623 x = points[i]; 624 pm2 = 0; 625 pm1 = 1; 626 pd2 = 0; 627 pd1 = 0; 628 pdd2 = 0; 629 pdd1 = 0; 630 k = 0; 631 if (degrees[k] == 0) { 632 if (B) B[i*ndegree+k] = pm1; 633 if (D) D[i*ndegree+k] = pd1; 634 if (D2) D2[i*ndegree+k] = pdd1; 635 k++; 636 } 637 for (j=1; j<=maxdegree; j++,k++) { 638 PetscReal p,d,dd; 639 p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 640 d = pd2 + (2*j-1)*pm1; 641 dd = pdd2 + (2*j-1)*pd1; 642 pm2 = pm1; 643 pm1 = p; 644 pd2 = pd1; 645 pd1 = d; 646 pdd2 = pdd1; 647 pdd1 = dd; 648 if (degrees[k] == j) { 649 if (B) B[i*ndegree+k] = p; 650 if (D) D[i*ndegree+k] = d; 651 if (D2) D2[i*ndegree+k] = dd; 652 } 653 } 654 } 655 PetscFunctionReturn(0); 656 } 657 658 /* if we can cobble together an eigenvalue / eigenvector routine for a symmetric tridiagonal system, we should use 659 * Golub & Welsch (Gauss-Jacobi) or Golub (Gauss-Lobatto-Jacobi) to compute Gaussian quadrature */ 660 #if (!defined(PETSC_MISSING_LAPACK_STEQR) || !defined(PETSC_MISSING_LAPACK_STEGR)) 661 #define PETSCDTGAUSSIANQUADRATURE_EIG 1 662 #endif 663 664 /* solve the symmetric tridiagonal eigenvalue system, writing the eigenvalues into eigs and the eigenvectors into V 665 * with lds n; diag and subdiag are overwritten */ 666 static PetscErrorCode PetscDTSymmetricTridiagonalEigensolve(PetscInt n, PetscReal diag[], PetscReal subdiag[], 667 PetscReal eigs[], PetscScalar V[]) 668 { 669 char jobz = 'V'; /* eigenvalues and eigenvectors */ 670 char range = 'A'; /* all eigenvalues will be found */ 671 PetscReal VL = 0.; /* ignored because range is 'A' */ 672 PetscReal VU = 0.; /* ignored because range is 'A' */ 673 PetscBLASInt IL = 0; /* ignored because range is 'A' */ 674 PetscBLASInt IU = 0; /* ignored because range is 'A' */ 675 PetscReal abstol = 0.; /* unused */ 676 PetscBLASInt bn, bm, ldz; /* bm will equal bn on exit */ 677 PetscBLASInt *isuppz; 678 PetscBLASInt lwork, liwork; 679 PetscReal workquery; 680 PetscBLASInt iworkquery; 681 PetscBLASInt *iwork; 682 PetscBLASInt info; 683 PetscReal *work = NULL; 684 PetscErrorCode ierr; 685 686 PetscFunctionBegin; 687 #if !defined(PETSCDTGAUSSIANQUADRATURE_EIG) 688 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found"); 689 #endif 690 ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr); 691 ierr = PetscBLASIntCast(n, &ldz);CHKERRQ(ierr); 692 #if !defined(PETSC_MISSING_LAPACK_STEGR) 693 ierr = PetscMalloc1(2 * n, &isuppz);CHKERRQ(ierr); 694 lwork = -1; 695 liwork = -1; 696 PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,&workquery,&lwork,&iworkquery,&liwork,&info)); 697 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error"); 698 lwork = (PetscBLASInt) workquery; 699 liwork = (PetscBLASInt) iworkquery; 700 ierr = PetscMalloc2(lwork, &work, liwork, &iwork);CHKERRQ(ierr); 701 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 702 PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,work,&lwork,iwork,&liwork,&info)); 703 ierr = PetscFPTrapPop();CHKERRQ(ierr); 704 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error"); 705 ierr = PetscFree2(work, iwork);CHKERRQ(ierr); 706 ierr = PetscFree(isuppz);CHKERRQ(ierr); 707 #elif !defined(PETSC_MISSING_LAPACK_STEQR) 708 jobz = 'I'; /* Compute eigenvalues and eigenvectors of the 709 tridiagonal matrix. Z is initialized to the identity 710 matrix. */ 711 ierr = PetscMalloc1(PetscMax(1,2*n-2),&work);CHKERRQ(ierr); 712 PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&bn,diag,subdiag,V,&ldz,work,&info)); 713 ierr = PetscFPTrapPop();CHKERRQ(ierr); 714 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 715 ierr = PetscFree(work);CHKERRQ(ierr); 716 ierr = PetscArraycpy(eigs,diag,n);CHKERRQ(ierr); 717 #endif 718 PetscFunctionReturn(0); 719 } 720 721 /* Formula for the weights at the endpoints (-1 and 1) of Gauss-Lobatto-Jacobi 722 * quadrature rules on the interval [-1, 1] */ 723 static PetscErrorCode PetscDTGaussLobattoJacobiEndweights_Internal(PetscInt n, PetscReal alpha, PetscReal beta, PetscReal *leftw, PetscReal *rightw) 724 { 725 PetscReal twoab1; 726 PetscInt m = n - 2; 727 PetscReal a = alpha + 1.; 728 PetscReal b = beta + 1.; 729 PetscReal gra, grb; 730 PetscErrorCode ierr; 731 732 PetscFunctionBegin; 733 twoab1 = PetscPowReal(2., a + b - 1.); 734 #if defined(PETSC_HAVE_LGAMMA) 735 grb = PetscExpReal(2. * PetscLGamma(b+1.) + PetscLGamma(m+1.) + PetscLGamma(m+a+1.) - 736 (PetscLGamma(m+b+1) + PetscLGamma(m+a+b+1.))); 737 gra = PetscExpReal(2. * PetscLGamma(a+1.) + PetscLGamma(m+1.) + PetscLGamma(m+b+1.) - 738 (PetscLGamma(m+a+1) + PetscLGamma(m+a+b+1.))); 739 #else 740 { 741 PetscInt alphai = (PetscInt) alpha; 742 PetscInt betai = (PetscInt) beta; 743 744 if ((PetscReal) alphai == alpha && (PetscReal) betai == beta) { 745 PetscReal binom1, binom2; 746 747 ierr = PetscDTBinomial(m+b, b, &binom1);CHKERRQ(ierr); 748 ierr = PetscDTBinomial(m+a+b, b, &binom2);CHKERRQ(ierr); 749 grb = 1./ (binom1 * binom2); 750 ierr = PetscDTBinomial(m+a, a, &binom1);CHKERRQ(ierr); 751 ierr = PetscDTBinomial(m+a+b, a, &binom2);CHKERRQ(ierr); 752 gra = 1./ (binom1 * binom2); 753 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable."); 754 } 755 #endif 756 *leftw = twoab1 * grb / b; 757 *rightw = twoab1 * gra / a; 758 PetscFunctionReturn(0); 759 } 760 761 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 762 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 763 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 764 { 765 PetscReal apb, pn1, pn2; 766 PetscInt k; 767 768 PetscFunctionBegin; 769 if (!n) {*P = 1.0; PetscFunctionReturn(0);} 770 if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 771 apb = a + b; 772 pn2 = 1.0; 773 pn1 = 0.5 * (a - b + (apb + 2.0) * x); 774 *P = 0.0; 775 for (k = 2; k < n+1; ++k) { 776 PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 777 PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 778 PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 779 PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 780 781 a2 = a2 / a1; 782 a3 = a3 / a1; 783 a4 = a4 / a1; 784 *P = (a2 + a3 * x) * pn1 - a4 * pn2; 785 pn2 = pn1; 786 pn1 = *P; 787 } 788 PetscFunctionReturn(0); 789 } 790 791 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 792 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscInt k, PetscReal *P) 793 { 794 PetscReal nP; 795 PetscInt i; 796 PetscErrorCode ierr; 797 798 PetscFunctionBegin; 799 if (k > n) {*P = 0.0; PetscFunctionReturn(0);} 800 ierr = PetscDTComputeJacobi(a+k, b+k, n-k, x, &nP);CHKERRQ(ierr); 801 for (i = 0; i < k; i++) nP *= (a + b + n + 1. + i) * 0.5; 802 *P = nP; 803 PetscFunctionReturn(0); 804 } 805 806 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 807 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 808 { 809 PetscFunctionBegin; 810 *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 811 *eta = y; 812 PetscFunctionReturn(0); 813 } 814 815 static PetscErrorCode PetscDTGaussJacobiQuadrature_Newton_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[]) 816 { 817 PetscInt maxIter = 100; 818 PetscReal eps = 1.0e-8; 819 PetscReal a1, a2, a3, a4, a5, a6; 820 PetscInt k; 821 PetscErrorCode ierr; 822 823 PetscFunctionBegin; 824 825 a1 = PetscPowReal(2.0, a+b+1); 826 #if defined(PETSC_HAVE_TGAMMA) 827 a2 = PetscTGamma(a + npoints + 1); 828 a3 = PetscTGamma(b + npoints + 1); 829 a4 = PetscTGamma(a + b + npoints + 1); 830 #else 831 { 832 PetscInt ia, ib; 833 834 ia = (PetscInt) a; 835 ib = (PetscInt) b; 836 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 */ 837 ierr = PetscDTFactorial(ia + npoints, &a2);CHKERRQ(ierr); 838 ierr = PetscDTFactorial(ib + npoints, &a3);CHKERRQ(ierr); 839 ierr = PetscDTFactorial(ia + ib + npoints, &a4);CHKERRQ(ierr); 840 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 841 } 842 #endif 843 844 ierr = PetscDTFactorial(npoints, &a5);CHKERRQ(ierr); 845 a6 = a1 * a2 * a3 / a4 / a5; 846 /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 847 Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 848 for (k = 0; k < npoints; ++k) { 849 PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 850 PetscInt j; 851 852 if (k > 0) r = 0.5 * (r + x[k-1]); 853 for (j = 0; j < maxIter; ++j) { 854 PetscReal s = 0.0, delta, f, fp; 855 PetscInt i; 856 857 for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 858 ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 859 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, 1, &fp);CHKERRQ(ierr); 860 delta = f / (fp - f * s); 861 r = r - delta; 862 if (PetscAbsReal(delta) < eps) break; 863 } 864 x[k] = r; 865 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], 1, &dP);CHKERRQ(ierr); 866 w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 867 } 868 PetscFunctionReturn(0); 869 } 870 871 /* Compute the diagonals of the Jacobi matrix used in Golub (& Welsch) algorithms for Gauss-(Lobatto)-Jacobi 872 * quadrature weight calculations on [-1,1] for exponents (1. + x)^a (1.-x)^b */ 873 static PetscErrorCode PetscDTJacobiMatrix_Internal(PetscInt nPoints, PetscReal a, PetscReal b, PetscReal *d, PetscReal *s) 874 { 875 PetscInt i; 876 877 PetscFunctionBegin; 878 for (i = 0; i < nPoints; i++) { 879 PetscReal A, B, Ap, C; 880 PetscInt si = i+1; 881 882 /* Jacobi three term recurrence */ 883 if (si > 1) { 884 A = (2.*si+a+b)*(2*si+a+b-1.) / (2.*si*(si+a+b)); 885 B = (a*a-b*b)*(2*si+a+b-1.) / (2.*si*(si+a+b)*(2*si+a+b-2)); 886 } else { 887 A = (a+b+2.) / 2.; 888 B = (a-b) / 2.; 889 } 890 Ap = (2*(si+1)+a+b)*(2*(si+1)+a+b-1) / (2*(si+1)*(si+1+a+b)); 891 C = (2*((si+1)+a-1)*((si+1)+b-1)*(2*(si+1)+a+b) / (2*(si+1)*((si+1)+a+b)*(2*(si+1)+a+b-2))); 892 d[i] = -B / A; 893 s[i] = (C / (A * Ap)); 894 } 895 PetscFunctionReturn(0); 896 } 897 898 static PetscErrorCode PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 899 { 900 PetscReal mu0; 901 PetscReal ga, gb, gab; 902 PetscInt i; 903 PetscErrorCode ierr; 904 905 PetscFunctionBegin; 906 ierr = PetscCitationsRegister(GolubWelschCitation, &GolubWelschCite);CHKERRQ(ierr); 907 908 #if defined(PETSC_HAVE_TGAMMA) 909 ga = PetscTGamma(a + 1); 910 gb = PetscTGamma(b + 1); 911 gab = PetscTGamma(a + b + 2); 912 #else 913 { 914 PetscInt ia, ib; 915 916 ia = (PetscInt) a; 917 ib = (PetscInt) b; 918 if (ia == a && ib == b && ia + 1 > 0 && ib + 1 > 0 && ia + ib + 2 > 0) { /* All gamma(x) terms are (x-1)! terms */ 919 ierr = PetscDTFactorial(ia, &ga);CHKERRQ(ierr); 920 ierr = PetscDTFactorial(ib, &gb);CHKERRQ(ierr); 921 ierr = PetscDTFactorial(ia + ib + 1, &gb);CHKERRQ(ierr); 922 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 923 } 924 #endif 925 mu0 = PetscPowReal(2.,a + b + 1.) * ga * gb / gab; 926 927 #if defined(PETSCDTGAUSSIANQUADRATURE_EIG) 928 { 929 PetscReal *diag, *subdiag; 930 PetscScalar *V; 931 932 ierr = PetscMalloc2(npoints, &diag, npoints, &subdiag);CHKERRQ(ierr); 933 ierr = PetscMalloc1(npoints*npoints, &V);CHKERRQ(ierr); 934 ierr = PetscDTJacobiMatrix_Internal(npoints, a, b, diag, subdiag);CHKERRQ(ierr); 935 for (i = 0; i < npoints - 1; i++) subdiag[i] = PetscSqrtReal(subdiag[i]); 936 ierr = PetscDTSymmetricTridiagonalEigensolve(npoints, diag, subdiag, x, V);CHKERRQ(ierr); 937 for (i = 0; i < npoints; i++) w[i] = PetscSqr(V[i * npoints]) * mu0; 938 ierr = PetscFree(V);CHKERRQ(ierr); 939 ierr = PetscFree2(diag, subdiag);CHKERRQ(ierr); 940 } 941 #else 942 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found"); 943 #endif 944 PetscFunctionReturn(0); 945 } 946 947 static PetscErrorCode PetscDTGaussJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton) 948 { 949 PetscErrorCode ierr; 950 951 PetscFunctionBegin; 952 if (npoints < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive"); 953 /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */ 954 if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1."); 955 if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1."); 956 957 if (newton) { 958 ierr = PetscDTGaussJacobiQuadrature_Newton_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr); 959 } else { 960 ierr = PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr); 961 } 962 if (alpha == beta) { /* symmetrize */ 963 PetscInt i; 964 for (i = 0; i < (npoints + 1) / 2; i++) { 965 PetscInt j = npoints - 1 - i; 966 PetscReal xi = x[i]; 967 PetscReal xj = x[j]; 968 PetscReal wi = w[i]; 969 PetscReal wj = w[j]; 970 971 x[i] = (xi - xj) / 2.; 972 x[j] = (xj - xi) / 2.; 973 w[i] = w[j] = (wi + wj) / 2.; 974 } 975 } 976 PetscFunctionReturn(0); 977 } 978 979 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[]) 980 { 981 PetscErrorCode ierr; 982 983 PetscFunctionBegin; 984 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, alpha, beta, x, w, gaussQuadratureNewton);CHKERRQ(ierr); 985 PetscFunctionReturn(0); 986 } 987 988 static PetscErrorCode PetscDTGaussLobattoJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton) 989 { 990 PetscInt i; 991 PetscErrorCode ierr; 992 993 PetscFunctionBegin; 994 if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive"); 995 /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */ 996 if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1."); 997 if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1."); 998 999 x[0] = -1.; 1000 x[npoints-1] = 1.; 1001 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, alpha+1., beta+1., &x[1], &w[1], newton);CHKERRQ(ierr); 1002 for (i = 1; i < npoints - 1; i++) { 1003 w[i] /= (1. - x[i]*x[i]); 1004 } 1005 ierr = PetscDTGaussLobattoJacobiEndweights_Internal(npoints, alpha, beta, &w[0], &w[npoints-1]);CHKERRQ(ierr); 1006 PetscFunctionReturn(0); 1007 } 1008 1009 /*@ 1010 PetscDTGaussQuadrature - create Gauss-Legendre quadrature 1011 1012 Not Collective 1013 1014 Input Arguments: 1015 + npoints - number of points 1016 . a - left end of interval (often-1) 1017 - b - right end of interval (often +1) 1018 1019 Output Arguments: 1020 + x - quadrature points 1021 - w - quadrature weights 1022 1023 Level: intermediate 1024 1025 References: 1026 . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 1027 1028 .seealso: PetscDTLegendreEval() 1029 @*/ 1030 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 1031 { 1032 PetscInt i; 1033 PetscErrorCode ierr; 1034 1035 PetscFunctionBegin; 1036 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, 0., 0., x, w, gaussQuadratureNewton);CHKERRQ(ierr); 1037 if (a != -1. || b != 1.) { 1038 for (i = 0; i < npoints; i++) { 1039 x[i] = (x[i] + 1.) * ((b - a) / 2.) + a; 1040 w[i] *= (b - a) / 2.; 1041 } 1042 } 1043 PetscFunctionReturn(0); 1044 } 1045 1046 /*@C 1047 PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre 1048 nodes of a given size on the domain [-1,1] 1049 1050 Not Collective 1051 1052 Input Parameter: 1053 + n - number of grid nodes 1054 - type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 1055 1056 Output Arguments: 1057 + x - quadrature points 1058 - w - quadrature weights 1059 1060 Notes: 1061 For n > 30 the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not 1062 close enough to the desired solution 1063 1064 These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes 1065 1066 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 1067 1068 Level: intermediate 1069 1070 .seealso: PetscDTGaussQuadrature() 1071 1072 @*/ 1073 PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w) 1074 { 1075 PetscBool newton; 1076 PetscErrorCode ierr; 1077 1078 PetscFunctionBegin; 1079 if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element"); 1080 newton = (PetscBool) (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA); 1081 ierr = PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, 0., 0., x, w, newton);CHKERRQ(ierr); 1082 PetscFunctionReturn(0); 1083 } 1084 1085 /*@ 1086 PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 1087 1088 Not Collective 1089 1090 Input Arguments: 1091 + dim - The spatial dimension 1092 . Nc - The number of components 1093 . npoints - number of points in one dimension 1094 . a - left end of interval (often-1) 1095 - b - right end of interval (often +1) 1096 1097 Output Argument: 1098 . q - A PetscQuadrature object 1099 1100 Level: intermediate 1101 1102 .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 1103 @*/ 1104 PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 1105 { 1106 PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c; 1107 PetscReal *x, *w, *xw, *ww; 1108 PetscErrorCode ierr; 1109 1110 PetscFunctionBegin; 1111 ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 1112 ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr); 1113 /* Set up the Golub-Welsch system */ 1114 switch (dim) { 1115 case 0: 1116 ierr = PetscFree(x);CHKERRQ(ierr); 1117 ierr = PetscFree(w);CHKERRQ(ierr); 1118 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 1119 ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 1120 x[0] = 0.0; 1121 for (c = 0; c < Nc; ++c) w[c] = 1.0; 1122 break; 1123 case 1: 1124 ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr); 1125 ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr); 1126 for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i]; 1127 ierr = PetscFree(ww);CHKERRQ(ierr); 1128 break; 1129 case 2: 1130 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 1131 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 1132 for (i = 0; i < npoints; ++i) { 1133 for (j = 0; j < npoints; ++j) { 1134 x[(i*npoints+j)*dim+0] = xw[i]; 1135 x[(i*npoints+j)*dim+1] = xw[j]; 1136 for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j]; 1137 } 1138 } 1139 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 1140 break; 1141 case 3: 1142 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 1143 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 1144 for (i = 0; i < npoints; ++i) { 1145 for (j = 0; j < npoints; ++j) { 1146 for (k = 0; k < npoints; ++k) { 1147 x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 1148 x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 1149 x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 1150 for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k]; 1151 } 1152 } 1153 } 1154 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 1155 break; 1156 default: 1157 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 1158 } 1159 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1160 ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 1161 ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 1162 ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr); 1163 PetscFunctionReturn(0); 1164 } 1165 1166 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 1167 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 1168 { 1169 PetscFunctionBegin; 1170 *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 1171 *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 1172 *zeta = z; 1173 PetscFunctionReturn(0); 1174 } 1175 1176 1177 /*@ 1178 PetscDTStroudConicalQuadrature - create Stroud conical quadrature for a simplex 1179 1180 Not Collective 1181 1182 Input Arguments: 1183 + dim - The simplex dimension 1184 . Nc - The number of components 1185 . npoints - The number of points in one dimension 1186 . a - left end of interval (often-1) 1187 - b - right end of interval (often +1) 1188 1189 Output Argument: 1190 . q - A PetscQuadrature object 1191 1192 Level: intermediate 1193 1194 References: 1195 . 1. - Karniadakis and Sherwin. FIAT 1196 1197 Note: For dim == 1, this is Gauss-Legendre quadrature 1198 1199 .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 1200 @*/ 1201 PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 1202 { 1203 PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints; 1204 PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 1205 PetscInt i, j, k, c; PetscErrorCode ierr; 1206 1207 PetscFunctionBegin; 1208 if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 1209 ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr); 1210 ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr); 1211 switch (dim) { 1212 case 0: 1213 ierr = PetscFree(x);CHKERRQ(ierr); 1214 ierr = PetscFree(w);CHKERRQ(ierr); 1215 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 1216 ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 1217 x[0] = 0.0; 1218 for (c = 0; c < Nc; ++c) w[c] = 1.0; 1219 break; 1220 case 1: 1221 ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr); 1222 ierr = PetscDTGaussJacobiQuadrature(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr); 1223 for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i]; 1224 ierr = PetscFree(wx);CHKERRQ(ierr); 1225 break; 1226 case 2: 1227 ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr); 1228 ierr = PetscDTGaussJacobiQuadrature(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 1229 ierr = PetscDTGaussJacobiQuadrature(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 1230 for (i = 0; i < npoints; ++i) { 1231 for (j = 0; j < npoints; ++j) { 1232 ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 1233 for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j]; 1234 } 1235 } 1236 ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 1237 break; 1238 case 3: 1239 ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr); 1240 ierr = PetscDTGaussJacobiQuadrature(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 1241 ierr = PetscDTGaussJacobiQuadrature(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 1242 ierr = PetscDTGaussJacobiQuadrature(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 1243 for (i = 0; i < npoints; ++i) { 1244 for (j = 0; j < npoints; ++j) { 1245 for (k = 0; k < npoints; ++k) { 1246 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); 1247 for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k]; 1248 } 1249 } 1250 } 1251 ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 1252 break; 1253 default: 1254 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 1255 } 1256 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1257 ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 1258 ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 1259 ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr); 1260 PetscFunctionReturn(0); 1261 } 1262 1263 /*@ 1264 PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 1265 1266 Not Collective 1267 1268 Input Arguments: 1269 + dim - The cell dimension 1270 . level - The number of points in one dimension, 2^l 1271 . a - left end of interval (often-1) 1272 - b - right end of interval (often +1) 1273 1274 Output Argument: 1275 . q - A PetscQuadrature object 1276 1277 Level: intermediate 1278 1279 .seealso: PetscDTGaussTensorQuadrature() 1280 @*/ 1281 PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 1282 { 1283 const PetscInt p = 16; /* Digits of precision in the evaluation */ 1284 const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1285 const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1286 const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 1287 PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 1288 PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 1289 PetscReal *x, *w; 1290 PetscInt K, k, npoints; 1291 PetscErrorCode ierr; 1292 1293 PetscFunctionBegin; 1294 if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 1295 if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 1296 /* Find K such that the weights are < 32 digits of precision */ 1297 for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 1298 wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 1299 } 1300 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1301 ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 1302 npoints = 2*K-1; 1303 ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 1304 ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 1305 /* Center term */ 1306 x[0] = beta; 1307 w[0] = 0.5*alpha*PETSC_PI; 1308 for (k = 1; k < K; ++k) { 1309 wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1310 xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h)); 1311 x[2*k-1] = -alpha*xk+beta; 1312 w[2*k-1] = wk; 1313 x[2*k+0] = alpha*xk+beta; 1314 w[2*k+0] = wk; 1315 } 1316 ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr); 1317 PetscFunctionReturn(0); 1318 } 1319 1320 PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1321 { 1322 const PetscInt p = 16; /* Digits of precision in the evaluation */ 1323 const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1324 const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1325 PetscReal h = 1.0; /* Step size, length between x_k */ 1326 PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 1327 PetscReal osum = 0.0; /* Integral on last level */ 1328 PetscReal psum = 0.0; /* Integral on the level before the last level */ 1329 PetscReal sum; /* Integral on current level */ 1330 PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 1331 PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 1332 PetscReal wk; /* Quadrature weight at x_k */ 1333 PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 1334 PetscInt d; /* Digits of precision in the integral */ 1335 1336 PetscFunctionBegin; 1337 if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 1338 /* Center term */ 1339 func(beta, &lval); 1340 sum = 0.5*alpha*PETSC_PI*lval; 1341 /* */ 1342 do { 1343 PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 1344 PetscInt k = 1; 1345 1346 ++l; 1347 /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 1348 /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 1349 psum = osum; 1350 osum = sum; 1351 h *= 0.5; 1352 sum *= 0.5; 1353 do { 1354 wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1355 yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1356 lx = -alpha*(1.0 - yk)+beta; 1357 rx = alpha*(1.0 - yk)+beta; 1358 func(lx, &lval); 1359 func(rx, &rval); 1360 lterm = alpha*wk*lval; 1361 maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 1362 sum += lterm; 1363 rterm = alpha*wk*rval; 1364 maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 1365 sum += rterm; 1366 ++k; 1367 /* Only need to evaluate every other point on refined levels */ 1368 if (l != 1) ++k; 1369 } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 1370 1371 d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 1372 d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 1373 d3 = PetscLog10Real(maxTerm) - p; 1374 if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 1375 else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 1376 d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1377 } while (d < digits && l < 12); 1378 *sol = sum; 1379 1380 PetscFunctionReturn(0); 1381 } 1382 1383 #if defined(PETSC_HAVE_MPFR) 1384 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1385 { 1386 const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 1387 PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 1388 mpfr_t alpha; /* Half-width of the integration interval */ 1389 mpfr_t beta; /* Center of the integration interval */ 1390 mpfr_t h; /* Step size, length between x_k */ 1391 mpfr_t osum; /* Integral on last level */ 1392 mpfr_t psum; /* Integral on the level before the last level */ 1393 mpfr_t sum; /* Integral on current level */ 1394 mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 1395 mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 1396 mpfr_t wk; /* Quadrature weight at x_k */ 1397 PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 1398 PetscInt d; /* Digits of precision in the integral */ 1399 mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 1400 1401 PetscFunctionBegin; 1402 if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 1403 /* Create high precision storage */ 1404 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); 1405 /* Initialization */ 1406 mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 1407 mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 1408 mpfr_set_d(osum, 0.0, MPFR_RNDN); 1409 mpfr_set_d(psum, 0.0, MPFR_RNDN); 1410 mpfr_set_d(h, 1.0, MPFR_RNDN); 1411 mpfr_const_pi(pi2, MPFR_RNDN); 1412 mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 1413 /* Center term */ 1414 func(0.5*(b+a), &lval); 1415 mpfr_set(sum, pi2, MPFR_RNDN); 1416 mpfr_mul(sum, sum, alpha, MPFR_RNDN); 1417 mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 1418 /* */ 1419 do { 1420 PetscReal d1, d2, d3, d4; 1421 PetscInt k = 1; 1422 1423 ++l; 1424 mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 1425 /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 1426 /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 1427 mpfr_set(psum, osum, MPFR_RNDN); 1428 mpfr_set(osum, sum, MPFR_RNDN); 1429 mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 1430 mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 1431 do { 1432 mpfr_set_si(kh, k, MPFR_RNDN); 1433 mpfr_mul(kh, kh, h, MPFR_RNDN); 1434 /* Weight */ 1435 mpfr_set(wk, h, MPFR_RNDN); 1436 mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 1437 mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 1438 mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 1439 mpfr_cosh(tmp, msinh, MPFR_RNDN); 1440 mpfr_sqr(tmp, tmp, MPFR_RNDN); 1441 mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 1442 mpfr_div(wk, wk, tmp, MPFR_RNDN); 1443 /* Abscissa */ 1444 mpfr_set_d(yk, 1.0, MPFR_RNDZ); 1445 mpfr_cosh(tmp, msinh, MPFR_RNDN); 1446 mpfr_div(yk, yk, tmp, MPFR_RNDZ); 1447 mpfr_exp(tmp, msinh, MPFR_RNDN); 1448 mpfr_div(yk, yk, tmp, MPFR_RNDZ); 1449 /* Quadrature points */ 1450 mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 1451 mpfr_mul(lx, lx, alpha, MPFR_RNDU); 1452 mpfr_add(lx, lx, beta, MPFR_RNDU); 1453 mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 1454 mpfr_mul(rx, rx, alpha, MPFR_RNDD); 1455 mpfr_add(rx, rx, beta, MPFR_RNDD); 1456 /* Evaluation */ 1457 func(mpfr_get_d(lx, MPFR_RNDU), &lval); 1458 func(mpfr_get_d(rx, MPFR_RNDD), &rval); 1459 /* Update */ 1460 mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 1461 mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 1462 mpfr_add(sum, sum, tmp, MPFR_RNDN); 1463 mpfr_abs(tmp, tmp, MPFR_RNDN); 1464 mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 1465 mpfr_set(curTerm, tmp, MPFR_RNDN); 1466 mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 1467 mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 1468 mpfr_add(sum, sum, tmp, MPFR_RNDN); 1469 mpfr_abs(tmp, tmp, MPFR_RNDN); 1470 mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 1471 mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 1472 ++k; 1473 /* Only need to evaluate every other point on refined levels */ 1474 if (l != 1) ++k; 1475 mpfr_log10(tmp, wk, MPFR_RNDN); 1476 mpfr_abs(tmp, tmp, MPFR_RNDN); 1477 } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 1478 mpfr_sub(tmp, sum, osum, MPFR_RNDN); 1479 mpfr_abs(tmp, tmp, MPFR_RNDN); 1480 mpfr_log10(tmp, tmp, MPFR_RNDN); 1481 d1 = mpfr_get_d(tmp, MPFR_RNDN); 1482 mpfr_sub(tmp, sum, psum, MPFR_RNDN); 1483 mpfr_abs(tmp, tmp, MPFR_RNDN); 1484 mpfr_log10(tmp, tmp, MPFR_RNDN); 1485 d2 = mpfr_get_d(tmp, MPFR_RNDN); 1486 mpfr_log10(tmp, maxTerm, MPFR_RNDN); 1487 d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 1488 mpfr_log10(tmp, curTerm, MPFR_RNDN); 1489 d4 = mpfr_get_d(tmp, MPFR_RNDN); 1490 d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1491 } while (d < digits && l < 8); 1492 *sol = mpfr_get_d(sum, MPFR_RNDN); 1493 /* Cleanup */ 1494 mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 1495 PetscFunctionReturn(0); 1496 } 1497 #else 1498 1499 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1500 { 1501 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 1502 } 1503 #endif 1504 1505 /* Overwrites A. Can only handle full-rank problems with m>=n 1506 * A in column-major format 1507 * Ainv in row-major format 1508 * tau has length m 1509 * worksize must be >= max(1,n) 1510 */ 1511 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1512 { 1513 PetscErrorCode ierr; 1514 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1515 PetscScalar *A,*Ainv,*R,*Q,Alpha; 1516 1517 PetscFunctionBegin; 1518 #if defined(PETSC_USE_COMPLEX) 1519 { 1520 PetscInt i,j; 1521 ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 1522 for (j=0; j<n; j++) { 1523 for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 1524 } 1525 mstride = m; 1526 } 1527 #else 1528 A = A_in; 1529 Ainv = Ainv_out; 1530 #endif 1531 1532 ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 1533 ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 1534 ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 1535 ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 1536 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1537 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1538 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1539 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1540 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1541 1542 /* Extract an explicit representation of Q */ 1543 Q = Ainv; 1544 ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr); 1545 K = N; /* full rank */ 1546 PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1547 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1548 1549 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1550 Alpha = 1.0; 1551 ldb = lda; 1552 PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 1553 /* Ainv is Q, overwritten with inverse */ 1554 1555 #if defined(PETSC_USE_COMPLEX) 1556 { 1557 PetscInt i; 1558 for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 1559 ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 1560 } 1561 #endif 1562 PetscFunctionReturn(0); 1563 } 1564 1565 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 1566 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 1567 { 1568 PetscErrorCode ierr; 1569 PetscReal *Bv; 1570 PetscInt i,j; 1571 1572 PetscFunctionBegin; 1573 ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 1574 /* Point evaluation of L_p on all the source vertices */ 1575 ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 1576 /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 1577 for (i=0; i<ninterval; i++) { 1578 for (j=0; j<ndegree; j++) { 1579 if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1580 else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1581 } 1582 } 1583 ierr = PetscFree(Bv);CHKERRQ(ierr); 1584 PetscFunctionReturn(0); 1585 } 1586 1587 /*@ 1588 PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 1589 1590 Not Collective 1591 1592 Input Arguments: 1593 + degree - degree of reconstruction polynomial 1594 . nsource - number of source intervals 1595 . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 1596 . ntarget - number of target intervals 1597 - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 1598 1599 Output Arguments: 1600 . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 1601 1602 Level: advanced 1603 1604 .seealso: PetscDTLegendreEval() 1605 @*/ 1606 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 1607 { 1608 PetscErrorCode ierr; 1609 PetscInt i,j,k,*bdegrees,worksize; 1610 PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 1611 PetscScalar *tau,*work; 1612 1613 PetscFunctionBegin; 1614 PetscValidRealPointer(sourcex,3); 1615 PetscValidRealPointer(targetx,5); 1616 PetscValidRealPointer(R,6); 1617 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); 1618 #if defined(PETSC_USE_DEBUG) 1619 for (i=0; i<nsource; i++) { 1620 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]); 1621 } 1622 for (i=0; i<ntarget; i++) { 1623 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]); 1624 } 1625 #endif 1626 xmin = PetscMin(sourcex[0],targetx[0]); 1627 xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 1628 center = (xmin + xmax)/2; 1629 hscale = (xmax - xmin)/2; 1630 worksize = nsource; 1631 ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 1632 ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 1633 for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 1634 for (i=0; i<=degree; i++) bdegrees[i] = i+1; 1635 ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 1636 ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 1637 for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 1638 ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 1639 for (i=0; i<ntarget; i++) { 1640 PetscReal rowsum = 0; 1641 for (j=0; j<nsource; j++) { 1642 PetscReal sum = 0; 1643 for (k=0; k<degree+1; k++) { 1644 sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 1645 } 1646 R[i*nsource+j] = sum; 1647 rowsum += sum; 1648 } 1649 for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 1650 } 1651 ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 1652 ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 1653 PetscFunctionReturn(0); 1654 } 1655 1656 /*@C 1657 PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points 1658 1659 Not Collective 1660 1661 Input Parameter: 1662 + n - the number of GLL nodes 1663 . nodes - the GLL nodes 1664 . weights - the GLL weights 1665 . f - the function values at the nodes 1666 1667 Output Parameter: 1668 . in - the value of the integral 1669 1670 Level: beginner 1671 1672 .seealso: PetscDTGaussLobattoLegendreQuadrature() 1673 1674 @*/ 1675 PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in) 1676 { 1677 PetscInt i; 1678 1679 PetscFunctionBegin; 1680 *in = 0.; 1681 for (i=0; i<n; i++) { 1682 *in += f[i]*f[i]*weights[i]; 1683 } 1684 PetscFunctionReturn(0); 1685 } 1686 1687 /*@C 1688 PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element 1689 1690 Not Collective 1691 1692 Input Parameter: 1693 + n - the number of GLL nodes 1694 . nodes - the GLL nodes 1695 . weights - the GLL weights 1696 1697 Output Parameter: 1698 . A - the stiffness element 1699 1700 Level: beginner 1701 1702 Notes: 1703 Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy() 1704 1705 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) 1706 1707 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1708 1709 @*/ 1710 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1711 { 1712 PetscReal **A; 1713 PetscErrorCode ierr; 1714 const PetscReal *gllnodes = nodes; 1715 const PetscInt p = n-1; 1716 PetscReal z0,z1,z2 = -1,x,Lpj,Lpr; 1717 PetscInt i,j,nn,r; 1718 1719 PetscFunctionBegin; 1720 ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1721 ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1722 for (i=1; i<n; i++) A[i] = A[i-1]+n; 1723 1724 for (j=1; j<p; j++) { 1725 x = gllnodes[j]; 1726 z0 = 1.; 1727 z1 = x; 1728 for (nn=1; nn<p; nn++) { 1729 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1730 z0 = z1; 1731 z1 = z2; 1732 } 1733 Lpj=z2; 1734 for (r=1; r<p; r++) { 1735 if (r == j) { 1736 A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj); 1737 } else { 1738 x = gllnodes[r]; 1739 z0 = 1.; 1740 z1 = x; 1741 for (nn=1; nn<p; nn++) { 1742 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1743 z0 = z1; 1744 z1 = z2; 1745 } 1746 Lpr = z2; 1747 A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r])); 1748 } 1749 } 1750 } 1751 for (j=1; j<p+1; j++) { 1752 x = gllnodes[j]; 1753 z0 = 1.; 1754 z1 = x; 1755 for (nn=1; nn<p; nn++) { 1756 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1757 z0 = z1; 1758 z1 = z2; 1759 } 1760 Lpj = z2; 1761 A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j])); 1762 A[0][j] = A[j][0]; 1763 } 1764 for (j=0; j<p; j++) { 1765 x = gllnodes[j]; 1766 z0 = 1.; 1767 z1 = x; 1768 for (nn=1; nn<p; nn++) { 1769 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1770 z0 = z1; 1771 z1 = z2; 1772 } 1773 Lpj=z2; 1774 1775 A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j])); 1776 A[j][p] = A[p][j]; 1777 } 1778 A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.; 1779 A[p][p]=A[0][0]; 1780 *AA = A; 1781 PetscFunctionReturn(0); 1782 } 1783 1784 /*@C 1785 PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element 1786 1787 Not Collective 1788 1789 Input Parameter: 1790 + n - the number of GLL nodes 1791 . nodes - the GLL nodes 1792 . weights - the GLL weightss 1793 - A - the stiffness element 1794 1795 Level: beginner 1796 1797 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate() 1798 1799 @*/ 1800 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1801 { 1802 PetscErrorCode ierr; 1803 1804 PetscFunctionBegin; 1805 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1806 ierr = PetscFree(*AA);CHKERRQ(ierr); 1807 *AA = NULL; 1808 PetscFunctionReturn(0); 1809 } 1810 1811 /*@C 1812 PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element 1813 1814 Not Collective 1815 1816 Input Parameter: 1817 + n - the number of GLL nodes 1818 . nodes - the GLL nodes 1819 . weights - the GLL weights 1820 1821 Output Parameter: 1822 . AA - the stiffness element 1823 - AAT - the transpose of AA (pass in NULL if you do not need this array) 1824 1825 Level: beginner 1826 1827 Notes: 1828 Destroy this with PetscGaussLobattoLegendreElementGradientDestroy() 1829 1830 You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1831 1832 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1833 1834 @*/ 1835 PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1836 { 1837 PetscReal **A, **AT = NULL; 1838 PetscErrorCode ierr; 1839 const PetscReal *gllnodes = nodes; 1840 const PetscInt p = n-1; 1841 PetscReal Li, Lj,d0; 1842 PetscInt i,j; 1843 1844 PetscFunctionBegin; 1845 ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1846 ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1847 for (i=1; i<n; i++) A[i] = A[i-1]+n; 1848 1849 if (AAT) { 1850 ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr); 1851 ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr); 1852 for (i=1; i<n; i++) AT[i] = AT[i-1]+n; 1853 } 1854 1855 if (n==1) {A[0][0] = 0.;} 1856 d0 = (PetscReal)p*((PetscReal)p+1.)/4.; 1857 for (i=0; i<n; i++) { 1858 for (j=0; j<n; j++) { 1859 A[i][j] = 0.; 1860 ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[i], &Li);CHKERRQ(ierr); 1861 ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[j], &Lj);CHKERRQ(ierr); 1862 if (i!=j) A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j])); 1863 if ((j==i) && (i==0)) A[i][j] = -d0; 1864 if (j==i && i==p) A[i][j] = d0; 1865 if (AT) AT[j][i] = A[i][j]; 1866 } 1867 } 1868 if (AAT) *AAT = AT; 1869 *AA = A; 1870 PetscFunctionReturn(0); 1871 } 1872 1873 /*@C 1874 PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate() 1875 1876 Not Collective 1877 1878 Input Parameter: 1879 + n - the number of GLL nodes 1880 . nodes - the GLL nodes 1881 . weights - the GLL weights 1882 . AA - the stiffness element 1883 - AAT - the transpose of the element 1884 1885 Level: beginner 1886 1887 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate() 1888 1889 @*/ 1890 PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1891 { 1892 PetscErrorCode ierr; 1893 1894 PetscFunctionBegin; 1895 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1896 ierr = PetscFree(*AA);CHKERRQ(ierr); 1897 *AA = NULL; 1898 if (*AAT) { 1899 ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr); 1900 ierr = PetscFree(*AAT);CHKERRQ(ierr); 1901 *AAT = NULL; 1902 } 1903 PetscFunctionReturn(0); 1904 } 1905 1906 /*@C 1907 PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element 1908 1909 Not Collective 1910 1911 Input Parameter: 1912 + n - the number of GLL nodes 1913 . nodes - the GLL nodes 1914 . weights - the GLL weightss 1915 1916 Output Parameter: 1917 . AA - the stiffness element 1918 1919 Level: beginner 1920 1921 Notes: 1922 Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy() 1923 1924 This is the same as the Gradient operator multiplied by the diagonal mass matrix 1925 1926 You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1927 1928 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy() 1929 1930 @*/ 1931 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1932 { 1933 PetscReal **D; 1934 PetscErrorCode ierr; 1935 const PetscReal *gllweights = weights; 1936 const PetscInt glln = n; 1937 PetscInt i,j; 1938 1939 PetscFunctionBegin; 1940 ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr); 1941 for (i=0; i<glln; i++){ 1942 for (j=0; j<glln; j++) { 1943 D[i][j] = gllweights[i]*D[i][j]; 1944 } 1945 } 1946 *AA = D; 1947 PetscFunctionReturn(0); 1948 } 1949 1950 /*@C 1951 PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element 1952 1953 Not Collective 1954 1955 Input Parameter: 1956 + n - the number of GLL nodes 1957 . nodes - the GLL nodes 1958 . weights - the GLL weights 1959 - A - advection 1960 1961 Level: beginner 1962 1963 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate() 1964 1965 @*/ 1966 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1967 { 1968 PetscErrorCode ierr; 1969 1970 PetscFunctionBegin; 1971 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1972 ierr = PetscFree(*AA);CHKERRQ(ierr); 1973 *AA = NULL; 1974 PetscFunctionReturn(0); 1975 } 1976 1977 PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1978 { 1979 PetscReal **A; 1980 PetscErrorCode ierr; 1981 const PetscReal *gllweights = weights; 1982 const PetscInt glln = n; 1983 PetscInt i,j; 1984 1985 PetscFunctionBegin; 1986 ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr); 1987 ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr); 1988 for (i=1; i<glln; i++) A[i] = A[i-1]+glln; 1989 if (glln==1) {A[0][0] = 0.;} 1990 for (i=0; i<glln; i++) { 1991 for (j=0; j<glln; j++) { 1992 A[i][j] = 0.; 1993 if (j==i) A[i][j] = gllweights[i]; 1994 } 1995 } 1996 *AA = A; 1997 PetscFunctionReturn(0); 1998 } 1999 2000 PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2001 { 2002 PetscErrorCode ierr; 2003 2004 PetscFunctionBegin; 2005 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 2006 ierr = PetscFree(*AA);CHKERRQ(ierr); 2007 *AA = NULL; 2008 PetscFunctionReturn(0); 2009 } 2010