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