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