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