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