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