1 /* Discretization tools */ 2 3 #include <petscconf.h> 4 #if defined(PETSC_HAVE_MATHIMF_H) 5 #include <mathimf.h> /* this needs to be included before math.h */ 6 #endif 7 #if defined(PETSC_HAVE_MPFR) 8 #include <mpfr.h> 9 #endif 10 11 #include <petscdt.h> /*I "petscdt.h" I*/ 12 #include <petscblaslapack.h> 13 #include <petsc/private/petscimpl.h> 14 #include <petsc/private/dtimpl.h> 15 #include <petscviewer.h> 16 #include <petscdmplex.h> 17 #include <petscdmshell.h> 18 19 static PetscBool GaussCite = PETSC_FALSE; 20 const char GaussCitation[] = "@article{GolubWelsch1969,\n" 21 " author = {Golub and Welsch},\n" 22 " title = {Calculation of Quadrature Rules},\n" 23 " journal = {Math. Comp.},\n" 24 " volume = {23},\n" 25 " number = {106},\n" 26 " pages = {221--230},\n" 27 " year = {1969}\n}\n"; 28 29 /*@ 30 PetscQuadratureCreate - Create a PetscQuadrature object 31 32 Collective on MPI_Comm 33 34 Input Parameter: 35 . comm - The communicator for the PetscQuadrature object 36 37 Output Parameter: 38 . q - The PetscQuadrature object 39 40 Level: beginner 41 42 .keywords: PetscQuadrature, quadrature, create 43 .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 44 @*/ 45 PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 46 { 47 PetscErrorCode ierr; 48 49 PetscFunctionBegin; 50 PetscValidPointer(q, 2); 51 ierr = PetscSysInitializePackage();CHKERRQ(ierr); 52 ierr = PetscHeaderCreate(*q,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 53 (*q)->dim = -1; 54 (*q)->Nc = 1; 55 (*q)->order = -1; 56 (*q)->numPoints = 0; 57 (*q)->points = NULL; 58 (*q)->weights = NULL; 59 PetscFunctionReturn(0); 60 } 61 62 /*@ 63 PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 64 65 Collective on PetscQuadrature 66 67 Input Parameter: 68 . q - The PetscQuadrature object 69 70 Output Parameter: 71 . r - The new PetscQuadrature object 72 73 Level: beginner 74 75 .keywords: PetscQuadrature, quadrature, clone 76 .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 77 @*/ 78 PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 79 { 80 PetscInt order, dim, Nc, Nq; 81 const PetscReal *points, *weights; 82 PetscReal *p, *w; 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 PetscValidPointer(q, 2); 87 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 88 ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 89 ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 90 ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 91 ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 92 ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr); 93 ierr = PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));CHKERRQ(ierr); 94 ierr = PetscMemcpy(w, weights, Nc * Nq * sizeof(PetscReal));CHKERRQ(ierr); 95 ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr); 96 PetscFunctionReturn(0); 97 } 98 99 /*@ 100 PetscQuadratureDestroy - Destroys a PetscQuadrature object 101 102 Collective on PetscQuadrature 103 104 Input Parameter: 105 . q - The PetscQuadrature object 106 107 Level: beginner 108 109 .keywords: PetscQuadrature, quadrature, destroy 110 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 111 @*/ 112 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 113 { 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 if (!*q) PetscFunctionReturn(0); 118 PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 119 if (--((PetscObject)(*q))->refct > 0) { 120 *q = NULL; 121 PetscFunctionReturn(0); 122 } 123 ierr = PetscFree((*q)->points);CHKERRQ(ierr); 124 ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 125 ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 126 PetscFunctionReturn(0); 127 } 128 129 /*@ 130 PetscQuadratureGetOrder - Return the order of the method 131 132 Not collective 133 134 Input Parameter: 135 . q - The PetscQuadrature object 136 137 Output Parameter: 138 . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 139 140 Level: intermediate 141 142 .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 143 @*/ 144 PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 145 { 146 PetscFunctionBegin; 147 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 148 PetscValidPointer(order, 2); 149 *order = q->order; 150 PetscFunctionReturn(0); 151 } 152 153 /*@ 154 PetscQuadratureSetOrder - Return the order of the method 155 156 Not collective 157 158 Input Parameters: 159 + q - The PetscQuadrature object 160 - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 161 162 Level: intermediate 163 164 .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 165 @*/ 166 PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 167 { 168 PetscFunctionBegin; 169 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 170 q->order = order; 171 PetscFunctionReturn(0); 172 } 173 174 /*@ 175 PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated 176 177 Not collective 178 179 Input Parameter: 180 . q - The PetscQuadrature object 181 182 Output Parameter: 183 . Nc - The number of components 184 185 Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 186 187 Level: intermediate 188 189 .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 190 @*/ 191 PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc) 192 { 193 PetscFunctionBegin; 194 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 195 PetscValidPointer(Nc, 2); 196 *Nc = q->Nc; 197 PetscFunctionReturn(0); 198 } 199 200 /*@ 201 PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated 202 203 Not collective 204 205 Input Parameters: 206 + q - The PetscQuadrature object 207 - Nc - The number of components 208 209 Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 210 211 Level: intermediate 212 213 .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 214 @*/ 215 PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc) 216 { 217 PetscFunctionBegin; 218 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 219 q->Nc = Nc; 220 PetscFunctionReturn(0); 221 } 222 223 /*@C 224 PetscQuadratureGetData - Returns the data defining the quadrature 225 226 Not collective 227 228 Input Parameter: 229 . q - The PetscQuadrature object 230 231 Output Parameters: 232 + dim - The spatial dimension 233 . Nc - The number of components 234 . npoints - The number of quadrature points 235 . points - The coordinates of each quadrature point 236 - weights - The weight of each quadrature point 237 238 Level: intermediate 239 240 Fortran Notes: 241 From Fortran you must call PetscQuadratureRestoreData() when you are done with the data 242 243 .keywords: PetscQuadrature, quadrature 244 .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 245 @*/ 246 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 247 { 248 PetscFunctionBegin; 249 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 250 if (dim) { 251 PetscValidPointer(dim, 2); 252 *dim = q->dim; 253 } 254 if (Nc) { 255 PetscValidPointer(Nc, 3); 256 *Nc = q->Nc; 257 } 258 if (npoints) { 259 PetscValidPointer(npoints, 4); 260 *npoints = q->numPoints; 261 } 262 if (points) { 263 PetscValidPointer(points, 5); 264 *points = q->points; 265 } 266 if (weights) { 267 PetscValidPointer(weights, 6); 268 *weights = q->weights; 269 } 270 PetscFunctionReturn(0); 271 } 272 273 /*@C 274 PetscQuadratureSetData - Sets the data defining the quadrature 275 276 Not collective 277 278 Input Parameters: 279 + q - The PetscQuadrature object 280 . dim - The spatial dimension 281 , Nc - The number of components 282 . npoints - The number of quadrature points 283 . points - The coordinates of each quadrature point 284 - weights - The weight of each quadrature point 285 286 Note: This routine owns the references to points and weights, so they msut be allocated using PetscMalloc() and the user should not free them. 287 288 Level: intermediate 289 290 .keywords: PetscQuadrature, quadrature 291 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 292 @*/ 293 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 294 { 295 PetscFunctionBegin; 296 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 297 if (dim >= 0) q->dim = dim; 298 if (Nc >= 0) q->Nc = Nc; 299 if (npoints >= 0) q->numPoints = npoints; 300 if (points) { 301 PetscValidPointer(points, 4); 302 q->points = points; 303 } 304 if (weights) { 305 PetscValidPointer(weights, 5); 306 q->weights = weights; 307 } 308 PetscFunctionReturn(0); 309 } 310 311 /*@C 312 PetscQuadratureView - Views a PetscQuadrature object 313 314 Collective on PetscQuadrature 315 316 Input Parameters: 317 + q - The PetscQuadrature object 318 - viewer - The PetscViewer object 319 320 Level: beginner 321 322 .keywords: PetscQuadrature, quadrature, view 323 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 324 @*/ 325 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 326 { 327 PetscInt q, d, c; 328 PetscErrorCode ierr; 329 330 PetscFunctionBegin; 331 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr); 332 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %D points with %D components\n (", quad->numPoints, quad->Nc);CHKERRQ(ierr);} 333 else {ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %D points\n (", quad->numPoints);CHKERRQ(ierr);} 334 for (q = 0; q < quad->numPoints; ++q) { 335 for (d = 0; d < quad->dim; ++d) { 336 if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 337 ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 338 } 339 if (quad->Nc > 1) { 340 ierr = PetscViewerASCIIPrintf(viewer, ") (");CHKERRQ(ierr); 341 for (c = 0; c < quad->Nc; ++c) { 342 if (c) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 343 ierr = PetscViewerASCIIPrintf(viewer, "%g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr); 344 } 345 ierr = PetscViewerASCIIPrintf(viewer, ")\n");CHKERRQ(ierr); 346 } else { 347 ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr); 348 } 349 } 350 PetscFunctionReturn(0); 351 } 352 353 /*@C 354 PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 355 356 Not collective 357 358 Input Parameter: 359 + q - The original PetscQuadrature 360 . numSubelements - The number of subelements the original element is divided into 361 . v0 - An array of the initial points for each subelement 362 - jac - An array of the Jacobian mappings from the reference to each subelement 363 364 Output Parameters: 365 . dim - The dimension 366 367 Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 368 369 Not available from Fortran 370 371 Level: intermediate 372 373 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 374 @*/ 375 PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 376 { 377 const PetscReal *points, *weights; 378 PetscReal *pointsRef, *weightsRef; 379 PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e; 380 PetscErrorCode ierr; 381 382 PetscFunctionBegin; 383 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 384 PetscValidPointer(v0, 3); 385 PetscValidPointer(jac, 4); 386 PetscValidPointer(qref, 5); 387 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 388 ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 389 ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 390 npointsRef = npoints*numSubelements; 391 ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 392 ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr); 393 for (c = 0; c < numSubelements; ++c) { 394 for (p = 0; p < npoints; ++p) { 395 for (d = 0; d < dim; ++d) { 396 pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 397 for (e = 0; e < dim; ++e) { 398 pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 399 } 400 } 401 /* Could also use detJ here */ 402 for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements; 403 } 404 } 405 ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 406 ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 407 PetscFunctionReturn(0); 408 } 409 410 /*@ 411 PetscDTLegendreEval - evaluate Legendre polynomial at points 412 413 Not Collective 414 415 Input Arguments: 416 + npoints - number of spatial points to evaluate at 417 . points - array of locations to evaluate at 418 . ndegree - number of basis degrees to evaluate 419 - degrees - sorted array of degrees to evaluate 420 421 Output Arguments: 422 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 423 . D - row-oriented derivative evaluation matrix (or NULL) 424 - D2 - row-oriented second derivative evaluation matrix (or NULL) 425 426 Level: intermediate 427 428 .seealso: PetscDTGaussQuadrature() 429 @*/ 430 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 431 { 432 PetscInt i,maxdegree; 433 434 PetscFunctionBegin; 435 if (!npoints || !ndegree) PetscFunctionReturn(0); 436 maxdegree = degrees[ndegree-1]; 437 for (i=0; i<npoints; i++) { 438 PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 439 PetscInt j,k; 440 x = points[i]; 441 pm2 = 0; 442 pm1 = 1; 443 pd2 = 0; 444 pd1 = 0; 445 pdd2 = 0; 446 pdd1 = 0; 447 k = 0; 448 if (degrees[k] == 0) { 449 if (B) B[i*ndegree+k] = pm1; 450 if (D) D[i*ndegree+k] = pd1; 451 if (D2) D2[i*ndegree+k] = pdd1; 452 k++; 453 } 454 for (j=1; j<=maxdegree; j++,k++) { 455 PetscReal p,d,dd; 456 p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 457 d = pd2 + (2*j-1)*pm1; 458 dd = pdd2 + (2*j-1)*pd1; 459 pm2 = pm1; 460 pm1 = p; 461 pd2 = pd1; 462 pd1 = d; 463 pdd2 = pdd1; 464 pdd1 = dd; 465 if (degrees[k] == j) { 466 if (B) B[i*ndegree+k] = p; 467 if (D) D[i*ndegree+k] = d; 468 if (D2) D2[i*ndegree+k] = dd; 469 } 470 } 471 } 472 PetscFunctionReturn(0); 473 } 474 475 /*@ 476 PetscDTGaussQuadrature - create Gauss quadrature 477 478 Not Collective 479 480 Input Arguments: 481 + npoints - number of points 482 . a - left end of interval (often-1) 483 - b - right end of interval (often +1) 484 485 Output Arguments: 486 + x - quadrature points 487 - w - quadrature weights 488 489 Level: intermediate 490 491 References: 492 . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 493 494 .seealso: PetscDTLegendreEval() 495 @*/ 496 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 497 { 498 PetscErrorCode ierr; 499 PetscInt i; 500 PetscReal *work; 501 PetscScalar *Z; 502 PetscBLASInt N,LDZ,info; 503 504 PetscFunctionBegin; 505 ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 506 /* Set up the Golub-Welsch system */ 507 for (i=0; i<npoints; i++) { 508 x[i] = 0; /* diagonal is 0 */ 509 if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 510 } 511 ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 512 ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 513 LDZ = N; 514 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 515 PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 516 ierr = PetscFPTrapPop();CHKERRQ(ierr); 517 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 518 519 for (i=0; i<(npoints+1)/2; i++) { 520 PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 521 x[i] = (a+b)/2 - y*(b-a)/2; 522 if (x[i] == -0.0) x[i] = 0.0; 523 x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 524 525 w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 526 } 527 ierr = PetscFree2(Z,work);CHKERRQ(ierr); 528 PetscFunctionReturn(0); 529 } 530 531 /*@ 532 PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 533 534 Not Collective 535 536 Input Arguments: 537 + dim - The spatial dimension 538 . Nc - The number of components 539 . npoints - number of points in one dimension 540 . a - left end of interval (often-1) 541 - b - right end of interval (often +1) 542 543 Output Argument: 544 . q - A PetscQuadrature object 545 546 Level: intermediate 547 548 .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 549 @*/ 550 PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 551 { 552 PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c; 553 PetscReal *x, *w, *xw, *ww; 554 PetscErrorCode ierr; 555 556 PetscFunctionBegin; 557 ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 558 ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr); 559 /* Set up the Golub-Welsch system */ 560 switch (dim) { 561 case 0: 562 ierr = PetscFree(x);CHKERRQ(ierr); 563 ierr = PetscFree(w);CHKERRQ(ierr); 564 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 565 ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 566 x[0] = 0.0; 567 for (c = 0; c < Nc; ++c) w[c] = 1.0; 568 break; 569 case 1: 570 ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr); 571 ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr); 572 for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i]; 573 ierr = PetscFree(ww);CHKERRQ(ierr); 574 break; 575 case 2: 576 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 577 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 578 for (i = 0; i < npoints; ++i) { 579 for (j = 0; j < npoints; ++j) { 580 x[(i*npoints+j)*dim+0] = xw[i]; 581 x[(i*npoints+j)*dim+1] = xw[j]; 582 for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j]; 583 } 584 } 585 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 586 break; 587 case 3: 588 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 589 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 590 for (i = 0; i < npoints; ++i) { 591 for (j = 0; j < npoints; ++j) { 592 for (k = 0; k < npoints; ++k) { 593 x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 594 x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 595 x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 596 for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k]; 597 } 598 } 599 } 600 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 601 break; 602 default: 603 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 604 } 605 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 606 ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 607 ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 608 PetscFunctionReturn(0); 609 } 610 611 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 612 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 613 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 614 { 615 PetscReal f = 1.0; 616 PetscInt i; 617 618 PetscFunctionBegin; 619 for (i = 1; i < n+1; ++i) f *= i; 620 *factorial = f; 621 PetscFunctionReturn(0); 622 } 623 624 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 625 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 626 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 627 { 628 PetscReal apb, pn1, pn2; 629 PetscInt k; 630 631 PetscFunctionBegin; 632 if (!n) {*P = 1.0; PetscFunctionReturn(0);} 633 if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 634 apb = a + b; 635 pn2 = 1.0; 636 pn1 = 0.5 * (a - b + (apb + 2.0) * x); 637 *P = 0.0; 638 for (k = 2; k < n+1; ++k) { 639 PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 640 PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 641 PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 642 PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 643 644 a2 = a2 / a1; 645 a3 = a3 / a1; 646 a4 = a4 / a1; 647 *P = (a2 + a3 * x) * pn1 - a4 * pn2; 648 pn2 = pn1; 649 pn1 = *P; 650 } 651 PetscFunctionReturn(0); 652 } 653 654 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 655 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 656 { 657 PetscReal nP; 658 PetscErrorCode ierr; 659 660 PetscFunctionBegin; 661 if (!n) {*P = 0.0; PetscFunctionReturn(0);} 662 ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 663 *P = 0.5 * (a + b + n + 1) * nP; 664 PetscFunctionReturn(0); 665 } 666 667 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 668 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 669 { 670 PetscFunctionBegin; 671 *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 672 *eta = y; 673 PetscFunctionReturn(0); 674 } 675 676 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 677 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 678 { 679 PetscFunctionBegin; 680 *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 681 *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 682 *zeta = z; 683 PetscFunctionReturn(0); 684 } 685 686 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 687 { 688 PetscInt maxIter = 100; 689 PetscReal eps = 1.0e-8; 690 PetscReal a1, a2, a3, a4, a5, a6; 691 PetscInt k; 692 PetscErrorCode ierr; 693 694 PetscFunctionBegin; 695 696 a1 = PetscPowReal(2.0, a+b+1); 697 #if defined(PETSC_HAVE_TGAMMA) 698 a2 = PetscTGamma(a + npoints + 1); 699 a3 = PetscTGamma(b + npoints + 1); 700 a4 = PetscTGamma(a + b + npoints + 1); 701 #else 702 { 703 PetscInt ia, ib; 704 705 ia = (PetscInt) a; 706 ib = (PetscInt) b; 707 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 */ 708 ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr); 709 ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr); 710 ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr); 711 } else { 712 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 713 } 714 } 715 #endif 716 717 ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 718 a6 = a1 * a2 * a3 / a4 / a5; 719 /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 720 Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 721 for (k = 0; k < npoints; ++k) { 722 PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 723 PetscInt j; 724 725 if (k > 0) r = 0.5 * (r + x[k-1]); 726 for (j = 0; j < maxIter; ++j) { 727 PetscReal s = 0.0, delta, f, fp; 728 PetscInt i; 729 730 for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 731 ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 732 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 733 delta = f / (fp - f * s); 734 r = r - delta; 735 if (PetscAbsReal(delta) < eps) break; 736 } 737 x[k] = r; 738 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 739 w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 740 } 741 PetscFunctionReturn(0); 742 } 743 744 /*@ 745 PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 746 747 Not Collective 748 749 Input Arguments: 750 + dim - The simplex dimension 751 . Nc - The number of components 752 . npoints - The number of points in one dimension 753 . a - left end of interval (often-1) 754 - b - right end of interval (often +1) 755 756 Output Argument: 757 . q - A PetscQuadrature object 758 759 Level: intermediate 760 761 References: 762 . 1. - Karniadakis and Sherwin. FIAT 763 764 .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 765 @*/ 766 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 767 { 768 PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints; 769 PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 770 PetscInt i, j, k, c; 771 PetscErrorCode ierr; 772 773 PetscFunctionBegin; 774 if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 775 ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr); 776 ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr); 777 switch (dim) { 778 case 0: 779 ierr = PetscFree(x);CHKERRQ(ierr); 780 ierr = PetscFree(w);CHKERRQ(ierr); 781 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 782 ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 783 x[0] = 0.0; 784 for (c = 0; c < Nc; ++c) w[c] = 1.0; 785 break; 786 case 1: 787 ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr); 788 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr); 789 for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i]; 790 ierr = PetscFree(wx);CHKERRQ(ierr); 791 break; 792 case 2: 793 ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr); 794 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 795 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 796 for (i = 0; i < npoints; ++i) { 797 for (j = 0; j < npoints; ++j) { 798 ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 799 for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j]; 800 } 801 } 802 ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 803 break; 804 case 3: 805 ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr); 806 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 807 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 808 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 809 for (i = 0; i < npoints; ++i) { 810 for (j = 0; j < npoints; ++j) { 811 for (k = 0; k < npoints; ++k) { 812 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); 813 for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k]; 814 } 815 } 816 } 817 ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 818 break; 819 default: 820 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 821 } 822 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 823 ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 824 ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 825 PetscFunctionReturn(0); 826 } 827 828 /*@ 829 PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 830 831 Not Collective 832 833 Input Arguments: 834 + dim - The cell dimension 835 . level - The number of points in one dimension, 2^l 836 . a - left end of interval (often-1) 837 - b - right end of interval (often +1) 838 839 Output Argument: 840 . q - A PetscQuadrature object 841 842 Level: intermediate 843 844 .seealso: PetscDTGaussTensorQuadrature() 845 @*/ 846 PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 847 { 848 const PetscInt p = 16; /* Digits of precision in the evaluation */ 849 const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 850 const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 851 const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 852 PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 853 PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 854 PetscReal *x, *w; 855 PetscInt K, k, npoints; 856 PetscErrorCode ierr; 857 858 PetscFunctionBegin; 859 if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 860 if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 861 /* Find K such that the weights are < 32 digits of precision */ 862 for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 863 wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 864 } 865 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 866 ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 867 npoints = 2*K-1; 868 ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 869 ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 870 /* Center term */ 871 x[0] = beta; 872 w[0] = 0.5*alpha*PETSC_PI; 873 for (k = 1; k < K; ++k) { 874 wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 875 xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h)); 876 x[2*k-1] = -alpha*xk+beta; 877 w[2*k-1] = wk; 878 x[2*k+0] = alpha*xk+beta; 879 w[2*k+0] = wk; 880 } 881 ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr); 882 PetscFunctionReturn(0); 883 } 884 885 PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 886 { 887 const PetscInt p = 16; /* Digits of precision in the evaluation */ 888 const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 889 const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 890 PetscReal h = 1.0; /* Step size, length between x_k */ 891 PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 892 PetscReal osum = 0.0; /* Integral on last level */ 893 PetscReal psum = 0.0; /* Integral on the level before the last level */ 894 PetscReal sum; /* Integral on current level */ 895 PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 896 PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 897 PetscReal wk; /* Quadrature weight at x_k */ 898 PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 899 PetscInt d; /* Digits of precision in the integral */ 900 901 PetscFunctionBegin; 902 if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 903 /* Center term */ 904 func(beta, &lval); 905 sum = 0.5*alpha*PETSC_PI*lval; 906 /* */ 907 do { 908 PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 909 PetscInt k = 1; 910 911 ++l; 912 /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 913 /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 914 psum = osum; 915 osum = sum; 916 h *= 0.5; 917 sum *= 0.5; 918 do { 919 wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 920 yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 921 lx = -alpha*(1.0 - yk)+beta; 922 rx = alpha*(1.0 - yk)+beta; 923 func(lx, &lval); 924 func(rx, &rval); 925 lterm = alpha*wk*lval; 926 maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 927 sum += lterm; 928 rterm = alpha*wk*rval; 929 maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 930 sum += rterm; 931 ++k; 932 /* Only need to evaluate every other point on refined levels */ 933 if (l != 1) ++k; 934 } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 935 936 d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 937 d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 938 d3 = PetscLog10Real(maxTerm) - p; 939 if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 940 else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 941 d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 942 } while (d < digits && l < 12); 943 *sol = sum; 944 945 PetscFunctionReturn(0); 946 } 947 948 #if defined(PETSC_HAVE_MPFR) 949 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 950 { 951 const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 952 PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 953 mpfr_t alpha; /* Half-width of the integration interval */ 954 mpfr_t beta; /* Center of the integration interval */ 955 mpfr_t h; /* Step size, length between x_k */ 956 mpfr_t osum; /* Integral on last level */ 957 mpfr_t psum; /* Integral on the level before the last level */ 958 mpfr_t sum; /* Integral on current level */ 959 mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 960 mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 961 mpfr_t wk; /* Quadrature weight at x_k */ 962 PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 963 PetscInt d; /* Digits of precision in the integral */ 964 mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 965 966 PetscFunctionBegin; 967 if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 968 /* Create high precision storage */ 969 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); 970 /* Initialization */ 971 mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 972 mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 973 mpfr_set_d(osum, 0.0, MPFR_RNDN); 974 mpfr_set_d(psum, 0.0, MPFR_RNDN); 975 mpfr_set_d(h, 1.0, MPFR_RNDN); 976 mpfr_const_pi(pi2, MPFR_RNDN); 977 mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 978 /* Center term */ 979 func(0.5*(b+a), &lval); 980 mpfr_set(sum, pi2, MPFR_RNDN); 981 mpfr_mul(sum, sum, alpha, MPFR_RNDN); 982 mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 983 /* */ 984 do { 985 PetscReal d1, d2, d3, d4; 986 PetscInt k = 1; 987 988 ++l; 989 mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 990 /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 991 /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 992 mpfr_set(psum, osum, MPFR_RNDN); 993 mpfr_set(osum, sum, MPFR_RNDN); 994 mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 995 mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 996 do { 997 mpfr_set_si(kh, k, MPFR_RNDN); 998 mpfr_mul(kh, kh, h, MPFR_RNDN); 999 /* Weight */ 1000 mpfr_set(wk, h, MPFR_RNDN); 1001 mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 1002 mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 1003 mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 1004 mpfr_cosh(tmp, msinh, MPFR_RNDN); 1005 mpfr_sqr(tmp, tmp, MPFR_RNDN); 1006 mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 1007 mpfr_div(wk, wk, tmp, MPFR_RNDN); 1008 /* Abscissa */ 1009 mpfr_set_d(yk, 1.0, MPFR_RNDZ); 1010 mpfr_cosh(tmp, msinh, MPFR_RNDN); 1011 mpfr_div(yk, yk, tmp, MPFR_RNDZ); 1012 mpfr_exp(tmp, msinh, MPFR_RNDN); 1013 mpfr_div(yk, yk, tmp, MPFR_RNDZ); 1014 /* Quadrature points */ 1015 mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 1016 mpfr_mul(lx, lx, alpha, MPFR_RNDU); 1017 mpfr_add(lx, lx, beta, MPFR_RNDU); 1018 mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 1019 mpfr_mul(rx, rx, alpha, MPFR_RNDD); 1020 mpfr_add(rx, rx, beta, MPFR_RNDD); 1021 /* Evaluation */ 1022 func(mpfr_get_d(lx, MPFR_RNDU), &lval); 1023 func(mpfr_get_d(rx, MPFR_RNDD), &rval); 1024 /* Update */ 1025 mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 1026 mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 1027 mpfr_add(sum, sum, tmp, MPFR_RNDN); 1028 mpfr_abs(tmp, tmp, MPFR_RNDN); 1029 mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 1030 mpfr_set(curTerm, tmp, MPFR_RNDN); 1031 mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 1032 mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 1033 mpfr_add(sum, sum, tmp, MPFR_RNDN); 1034 mpfr_abs(tmp, tmp, MPFR_RNDN); 1035 mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 1036 mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 1037 ++k; 1038 /* Only need to evaluate every other point on refined levels */ 1039 if (l != 1) ++k; 1040 mpfr_log10(tmp, wk, MPFR_RNDN); 1041 mpfr_abs(tmp, tmp, MPFR_RNDN); 1042 } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 1043 mpfr_sub(tmp, sum, osum, MPFR_RNDN); 1044 mpfr_abs(tmp, tmp, MPFR_RNDN); 1045 mpfr_log10(tmp, tmp, MPFR_RNDN); 1046 d1 = mpfr_get_d(tmp, MPFR_RNDN); 1047 mpfr_sub(tmp, sum, psum, MPFR_RNDN); 1048 mpfr_abs(tmp, tmp, MPFR_RNDN); 1049 mpfr_log10(tmp, tmp, MPFR_RNDN); 1050 d2 = mpfr_get_d(tmp, MPFR_RNDN); 1051 mpfr_log10(tmp, maxTerm, MPFR_RNDN); 1052 d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 1053 mpfr_log10(tmp, curTerm, MPFR_RNDN); 1054 d4 = mpfr_get_d(tmp, MPFR_RNDN); 1055 d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1056 } while (d < digits && l < 8); 1057 *sol = mpfr_get_d(sum, MPFR_RNDN); 1058 /* Cleanup */ 1059 mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 1060 PetscFunctionReturn(0); 1061 } 1062 #else 1063 1064 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1065 { 1066 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 1067 } 1068 #endif 1069 1070 /* Overwrites A. Can only handle full-rank problems with m>=n 1071 * A in column-major format 1072 * Ainv in row-major format 1073 * tau has length m 1074 * worksize must be >= max(1,n) 1075 */ 1076 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1077 { 1078 PetscErrorCode ierr; 1079 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1080 PetscScalar *A,*Ainv,*R,*Q,Alpha; 1081 1082 PetscFunctionBegin; 1083 #if defined(PETSC_USE_COMPLEX) 1084 { 1085 PetscInt i,j; 1086 ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 1087 for (j=0; j<n; j++) { 1088 for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 1089 } 1090 mstride = m; 1091 } 1092 #else 1093 A = A_in; 1094 Ainv = Ainv_out; 1095 #endif 1096 1097 ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 1098 ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 1099 ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 1100 ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 1101 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1102 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1103 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1104 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1105 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1106 1107 /* Extract an explicit representation of Q */ 1108 Q = Ainv; 1109 ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 1110 K = N; /* full rank */ 1111 PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1112 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1113 1114 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1115 Alpha = 1.0; 1116 ldb = lda; 1117 PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 1118 /* Ainv is Q, overwritten with inverse */ 1119 1120 #if defined(PETSC_USE_COMPLEX) 1121 { 1122 PetscInt i; 1123 for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 1124 ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 1125 } 1126 #endif 1127 PetscFunctionReturn(0); 1128 } 1129 1130 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 1131 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 1132 { 1133 PetscErrorCode ierr; 1134 PetscReal *Bv; 1135 PetscInt i,j; 1136 1137 PetscFunctionBegin; 1138 ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 1139 /* Point evaluation of L_p on all the source vertices */ 1140 ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 1141 /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 1142 for (i=0; i<ninterval; i++) { 1143 for (j=0; j<ndegree; j++) { 1144 if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1145 else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1146 } 1147 } 1148 ierr = PetscFree(Bv);CHKERRQ(ierr); 1149 PetscFunctionReturn(0); 1150 } 1151 1152 /*@ 1153 PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 1154 1155 Not Collective 1156 1157 Input Arguments: 1158 + degree - degree of reconstruction polynomial 1159 . nsource - number of source intervals 1160 . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 1161 . ntarget - number of target intervals 1162 - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 1163 1164 Output Arguments: 1165 . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 1166 1167 Level: advanced 1168 1169 .seealso: PetscDTLegendreEval() 1170 @*/ 1171 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 1172 { 1173 PetscErrorCode ierr; 1174 PetscInt i,j,k,*bdegrees,worksize; 1175 PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 1176 PetscScalar *tau,*work; 1177 1178 PetscFunctionBegin; 1179 PetscValidRealPointer(sourcex,3); 1180 PetscValidRealPointer(targetx,5); 1181 PetscValidRealPointer(R,6); 1182 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); 1183 #if defined(PETSC_USE_DEBUG) 1184 for (i=0; i<nsource; i++) { 1185 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]); 1186 } 1187 for (i=0; i<ntarget; i++) { 1188 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]); 1189 } 1190 #endif 1191 xmin = PetscMin(sourcex[0],targetx[0]); 1192 xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 1193 center = (xmin + xmax)/2; 1194 hscale = (xmax - xmin)/2; 1195 worksize = nsource; 1196 ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 1197 ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 1198 for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 1199 for (i=0; i<=degree; i++) bdegrees[i] = i+1; 1200 ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 1201 ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 1202 for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 1203 ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 1204 for (i=0; i<ntarget; i++) { 1205 PetscReal rowsum = 0; 1206 for (j=0; j<nsource; j++) { 1207 PetscReal sum = 0; 1208 for (k=0; k<degree+1; k++) { 1209 sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 1210 } 1211 R[i*nsource+j] = sum; 1212 rowsum += sum; 1213 } 1214 for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 1215 } 1216 ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 1217 ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 1218 PetscFunctionReturn(0); 1219 } 1220