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