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,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 52 (*q)->dim = -1; 53 (*q)->numPoints = 0; 54 (*q)->points = NULL; 55 (*q)->weights = NULL; 56 PetscFunctionReturn(0); 57 } 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "PetscQuadratureDestroy" 61 /*@ 62 PetscQuadratureDestroy - Destroys a PetscQuadrature object 63 64 Collective on PetscQuadrature 65 66 Input Parameter: 67 . q - The PetscQuadrature object 68 69 Level: beginner 70 71 .keywords: PetscQuadrature, quadrature, destroy 72 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 73 @*/ 74 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 75 { 76 PetscErrorCode ierr; 77 78 PetscFunctionBegin; 79 if (!*q) PetscFunctionReturn(0); 80 PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); 81 if (--((PetscObject)(*q))->refct > 0) { 82 *q = NULL; 83 PetscFunctionReturn(0); 84 } 85 ierr = PetscFree((*q)->points);CHKERRQ(ierr); 86 ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 87 ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 88 PetscFunctionReturn(0); 89 } 90 91 #undef __FUNCT__ 92 #define __FUNCT__ "PetscQuadratureGetData" 93 /*@C 94 PetscQuadratureGetData - Returns the data defining the quadrature 95 96 Not collective 97 98 Input Parameter: 99 . q - The PetscQuadrature object 100 101 Output Parameters: 102 + dim - The spatial dimension 103 . npoints - The number of quadrature points 104 . points - The coordinates of each quadrature point 105 - weights - The weight of each quadrature point 106 107 Level: intermediate 108 109 .keywords: PetscQuadrature, quadrature 110 .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 111 @*/ 112 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 113 { 114 PetscFunctionBegin; 115 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 116 if (dim) { 117 PetscValidPointer(dim, 2); 118 *dim = q->dim; 119 } 120 if (npoints) { 121 PetscValidPointer(npoints, 3); 122 *npoints = q->numPoints; 123 } 124 if (points) { 125 PetscValidPointer(points, 4); 126 *points = q->points; 127 } 128 if (weights) { 129 PetscValidPointer(weights, 5); 130 *weights = q->weights; 131 } 132 PetscFunctionReturn(0); 133 } 134 135 #undef __FUNCT__ 136 #define __FUNCT__ "PetscQuadratureSetData" 137 /*@C 138 PetscQuadratureSetData - Sets the data defining the quadrature 139 140 Not collective 141 142 Input Parameters: 143 + q - The PetscQuadrature object 144 . dim - The spatial dimension 145 . npoints - The number of quadrature points 146 . points - The coordinates of each quadrature point 147 - weights - The weight of each quadrature point 148 149 Level: intermediate 150 151 .keywords: PetscQuadrature, quadrature 152 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 153 @*/ 154 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 155 { 156 PetscFunctionBegin; 157 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 158 if (dim >= 0) q->dim = dim; 159 if (npoints >= 0) q->numPoints = npoints; 160 if (points) { 161 PetscValidPointer(points, 4); 162 q->points = points; 163 } 164 if (weights) { 165 PetscValidPointer(weights, 5); 166 q->weights = weights; 167 } 168 PetscFunctionReturn(0); 169 } 170 171 #undef __FUNCT__ 172 #define __FUNCT__ "PetscQuadratureView" 173 /*@C 174 PetscQuadratureView - Views a PetscQuadrature object 175 176 Collective on PetscQuadrature 177 178 Input Parameters: 179 + q - The PetscQuadrature object 180 - viewer - The PetscViewer object 181 182 Level: beginner 183 184 .keywords: PetscQuadrature, quadrature, view 185 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 186 @*/ 187 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 188 { 189 PetscInt q, d; 190 PetscErrorCode ierr; 191 192 PetscFunctionBegin; 193 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr); 194 ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n (", quad->numPoints);CHKERRQ(ierr); 195 for (q = 0; q < quad->numPoints; ++q) { 196 for (d = 0; d < quad->dim; ++d) { 197 if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); 198 ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 199 } 200 ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr); 201 } 202 PetscFunctionReturn(0); 203 } 204 205 #undef __FUNCT__ 206 #define __FUNCT__ "PetscDTLegendreEval" 207 /*@ 208 PetscDTLegendreEval - evaluate Legendre polynomial at points 209 210 Not Collective 211 212 Input Arguments: 213 + npoints - number of spatial points to evaluate at 214 . points - array of locations to evaluate at 215 . ndegree - number of basis degrees to evaluate 216 - degrees - sorted array of degrees to evaluate 217 218 Output Arguments: 219 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 220 . D - row-oriented derivative evaluation matrix (or NULL) 221 - D2 - row-oriented second derivative evaluation matrix (or NULL) 222 223 Level: intermediate 224 225 .seealso: PetscDTGaussQuadrature() 226 @*/ 227 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 228 { 229 PetscInt i,maxdegree; 230 231 PetscFunctionBegin; 232 if (!npoints || !ndegree) PetscFunctionReturn(0); 233 maxdegree = degrees[ndegree-1]; 234 for (i=0; i<npoints; i++) { 235 PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 236 PetscInt j,k; 237 x = points[i]; 238 pm2 = 0; 239 pm1 = 1; 240 pd2 = 0; 241 pd1 = 0; 242 pdd2 = 0; 243 pdd1 = 0; 244 k = 0; 245 if (degrees[k] == 0) { 246 if (B) B[i*ndegree+k] = pm1; 247 if (D) D[i*ndegree+k] = pd1; 248 if (D2) D2[i*ndegree+k] = pdd1; 249 k++; 250 } 251 for (j=1; j<=maxdegree; j++,k++) { 252 PetscReal p,d,dd; 253 p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 254 d = pd2 + (2*j-1)*pm1; 255 dd = pdd2 + (2*j-1)*pd1; 256 pm2 = pm1; 257 pm1 = p; 258 pd2 = pd1; 259 pd1 = d; 260 pdd2 = pdd1; 261 pdd1 = dd; 262 if (degrees[k] == j) { 263 if (B) B[i*ndegree+k] = p; 264 if (D) D[i*ndegree+k] = d; 265 if (D2) D2[i*ndegree+k] = dd; 266 } 267 } 268 } 269 PetscFunctionReturn(0); 270 } 271 272 #undef __FUNCT__ 273 #define __FUNCT__ "PetscDTGaussQuadrature" 274 /*@ 275 PetscDTGaussQuadrature - create Gauss quadrature 276 277 Not Collective 278 279 Input Arguments: 280 + npoints - number of points 281 . a - left end of interval (often-1) 282 - b - right end of interval (often +1) 283 284 Output Arguments: 285 + x - quadrature points 286 - w - quadrature weights 287 288 Level: intermediate 289 290 References: 291 Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 292 293 .seealso: PetscDTLegendreEval() 294 @*/ 295 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 296 { 297 PetscErrorCode ierr; 298 PetscInt i; 299 PetscReal *work; 300 PetscScalar *Z; 301 PetscBLASInt N,LDZ,info; 302 303 PetscFunctionBegin; 304 ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr); 305 /* Set up the Golub-Welsch system */ 306 for (i=0; i<npoints; i++) { 307 x[i] = 0; /* diagonal is 0 */ 308 if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 309 } 310 ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 311 ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 312 LDZ = N; 313 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 314 PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 315 ierr = PetscFPTrapPop();CHKERRQ(ierr); 316 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 317 318 for (i=0; i<(npoints+1)/2; i++) { 319 PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 320 x[i] = (a+b)/2 - y*(b-a)/2; 321 x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 322 323 w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints]))); 324 } 325 ierr = PetscFree2(Z,work);CHKERRQ(ierr); 326 PetscFunctionReturn(0); 327 } 328 329 #undef __FUNCT__ 330 #define __FUNCT__ "PetscDTGaussTensorQuadrature" 331 /*@ 332 PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 333 334 Not Collective 335 336 Input Arguments: 337 + dim - The spatial dimension 338 . npoints - number of points in one dimension 339 . a - left end of interval (often-1) 340 - b - right end of interval (often +1) 341 342 Output Argument: 343 . q - A PetscQuadrature object 344 345 Level: intermediate 346 347 .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 348 @*/ 349 PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 350 { 351 PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k; 352 PetscReal *x, *w, *xw, *ww; 353 PetscErrorCode ierr; 354 355 PetscFunctionBegin; 356 ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 357 ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr); 358 /* Set up the Golub-Welsch system */ 359 switch (dim) { 360 case 0: 361 ierr = PetscFree(x);CHKERRQ(ierr); 362 ierr = PetscFree(w);CHKERRQ(ierr); 363 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 364 ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 365 x[0] = 0.0; 366 w[0] = 1.0; 367 break; 368 case 1: 369 ierr = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr); 370 break; 371 case 2: 372 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 373 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 374 for (i = 0; i < npoints; ++i) { 375 for (j = 0; j < npoints; ++j) { 376 x[(i*npoints+j)*dim+0] = xw[i]; 377 x[(i*npoints+j)*dim+1] = xw[j]; 378 w[i*npoints+j] = ww[i] * ww[j]; 379 } 380 } 381 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 382 break; 383 case 3: 384 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 385 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 386 for (i = 0; i < npoints; ++i) { 387 for (j = 0; j < npoints; ++j) { 388 for (k = 0; k < npoints; ++k) { 389 x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 390 x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 391 x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 392 w[(i*npoints+j)*npoints+k] = ww[i] * ww[j] * ww[k]; 393 } 394 } 395 } 396 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 397 break; 398 default: 399 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 400 } 401 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 402 ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr); 403 PetscFunctionReturn(0); 404 } 405 406 #undef __FUNCT__ 407 #define __FUNCT__ "PetscDTFactorial_Internal" 408 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 409 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 410 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 411 { 412 PetscReal f = 1.0; 413 PetscInt i; 414 415 PetscFunctionBegin; 416 for (i = 1; i < n+1; ++i) f *= i; 417 *factorial = f; 418 PetscFunctionReturn(0); 419 } 420 421 #undef __FUNCT__ 422 #define __FUNCT__ "PetscDTComputeJacobi" 423 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 424 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 425 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 426 { 427 PetscReal apb, pn1, pn2; 428 PetscInt k; 429 430 PetscFunctionBegin; 431 if (!n) {*P = 1.0; PetscFunctionReturn(0);} 432 if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 433 apb = a + b; 434 pn2 = 1.0; 435 pn1 = 0.5 * (a - b + (apb + 2.0) * x); 436 *P = 0.0; 437 for (k = 2; k < n+1; ++k) { 438 PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 439 PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 440 PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 441 PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 442 443 a2 = a2 / a1; 444 a3 = a3 / a1; 445 a4 = a4 / a1; 446 *P = (a2 + a3 * x) * pn1 - a4 * pn2; 447 pn2 = pn1; 448 pn1 = *P; 449 } 450 PetscFunctionReturn(0); 451 } 452 453 #undef __FUNCT__ 454 #define __FUNCT__ "PetscDTComputeJacobiDerivative" 455 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 456 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 457 { 458 PetscReal nP; 459 PetscErrorCode ierr; 460 461 PetscFunctionBegin; 462 if (!n) {*P = 0.0; PetscFunctionReturn(0);} 463 ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 464 *P = 0.5 * (a + b + n + 1) * nP; 465 PetscFunctionReturn(0); 466 } 467 468 #undef __FUNCT__ 469 #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 470 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 471 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 472 { 473 PetscFunctionBegin; 474 *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 475 *eta = y; 476 PetscFunctionReturn(0); 477 } 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 481 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 482 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 483 { 484 PetscFunctionBegin; 485 *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 486 *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 487 *zeta = z; 488 PetscFunctionReturn(0); 489 } 490 491 #undef __FUNCT__ 492 #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 493 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 494 { 495 PetscInt maxIter = 100; 496 PetscReal eps = 1.0e-8; 497 PetscReal a1, a2, a3, a4, a5, a6; 498 PetscInt k; 499 PetscErrorCode ierr; 500 501 PetscFunctionBegin; 502 503 a1 = PetscPowReal(2.0, a+b+1); 504 #if defined(PETSC_HAVE_TGAMMA) 505 a2 = PetscTGamma(a + npoints + 1); 506 a3 = PetscTGamma(b + npoints + 1); 507 a4 = PetscTGamma(a + b + npoints + 1); 508 #else 509 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 510 #endif 511 512 ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 513 a6 = a1 * a2 * a3 / a4 / a5; 514 /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 515 Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 516 for (k = 0; k < npoints; ++k) { 517 PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 518 PetscInt j; 519 520 if (k > 0) r = 0.5 * (r + x[k-1]); 521 for (j = 0; j < maxIter; ++j) { 522 PetscReal s = 0.0, delta, f, fp; 523 PetscInt i; 524 525 for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 526 ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 527 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 528 delta = f / (fp - f * s); 529 r = r - delta; 530 if (PetscAbsReal(delta) < eps) break; 531 } 532 x[k] = r; 533 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 534 w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 535 } 536 PetscFunctionReturn(0); 537 } 538 539 #undef __FUNCT__ 540 #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 541 /*@C 542 PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 543 544 Not Collective 545 546 Input Arguments: 547 + dim - The simplex dimension 548 . order - The number of points in one dimension 549 . a - left end of interval (often-1) 550 - b - right end of interval (often +1) 551 552 Output Argument: 553 . q - A PetscQuadrature object 554 555 Level: intermediate 556 557 References: 558 Karniadakis and Sherwin. 559 FIAT 560 561 .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 562 @*/ 563 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q) 564 { 565 PetscInt npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order; 566 PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 567 PetscInt i, j, k; 568 PetscErrorCode ierr; 569 570 PetscFunctionBegin; 571 if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 572 ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 573 ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 574 switch (dim) { 575 case 0: 576 ierr = PetscFree(x);CHKERRQ(ierr); 577 ierr = PetscFree(w);CHKERRQ(ierr); 578 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 579 ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 580 x[0] = 0.0; 581 w[0] = 1.0; 582 break; 583 case 1: 584 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr); 585 break; 586 case 2: 587 ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr); 588 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 589 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 590 for (i = 0; i < order; ++i) { 591 for (j = 0; j < order; ++j) { 592 ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr); 593 w[i*order+j] = 0.5 * wx[i] * wy[j]; 594 } 595 } 596 ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 597 break; 598 case 3: 599 ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr); 600 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 601 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 602 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 603 for (i = 0; i < order; ++i) { 604 for (j = 0; j < order; ++j) { 605 for (k = 0; k < order; ++k) { 606 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); 607 w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k]; 608 } 609 } 610 } 611 ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 612 break; 613 default: 614 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 615 } 616 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 617 ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr); 618 PetscFunctionReturn(0); 619 } 620 621 #undef __FUNCT__ 622 #define __FUNCT__ "PetscDTPseudoInverseQR" 623 /* Overwrites A. Can only handle full-rank problems with m>=n 624 * A in column-major format 625 * Ainv in row-major format 626 * tau has length m 627 * worksize must be >= max(1,n) 628 */ 629 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 630 { 631 PetscErrorCode ierr; 632 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 633 PetscScalar *A,*Ainv,*R,*Q,Alpha; 634 635 PetscFunctionBegin; 636 #if defined(PETSC_USE_COMPLEX) 637 { 638 PetscInt i,j; 639 ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 640 for (j=0; j<n; j++) { 641 for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 642 } 643 mstride = m; 644 } 645 #else 646 A = A_in; 647 Ainv = Ainv_out; 648 #endif 649 650 ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 651 ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 652 ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 653 ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 654 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 655 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 656 ierr = PetscFPTrapPop();CHKERRQ(ierr); 657 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 658 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 659 660 /* Extract an explicit representation of Q */ 661 Q = Ainv; 662 ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 663 K = N; /* full rank */ 664 PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 665 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 666 667 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 668 Alpha = 1.0; 669 ldb = lda; 670 PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 671 /* Ainv is Q, overwritten with inverse */ 672 673 #if defined(PETSC_USE_COMPLEX) 674 { 675 PetscInt i; 676 for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 677 ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 678 } 679 #endif 680 PetscFunctionReturn(0); 681 } 682 683 #undef __FUNCT__ 684 #define __FUNCT__ "PetscDTLegendreIntegrate" 685 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 686 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 687 { 688 PetscErrorCode ierr; 689 PetscReal *Bv; 690 PetscInt i,j; 691 692 PetscFunctionBegin; 693 ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 694 /* Point evaluation of L_p on all the source vertices */ 695 ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 696 /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 697 for (i=0; i<ninterval; i++) { 698 for (j=0; j<ndegree; j++) { 699 if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 700 else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 701 } 702 } 703 ierr = PetscFree(Bv);CHKERRQ(ierr); 704 PetscFunctionReturn(0); 705 } 706 707 #undef __FUNCT__ 708 #define __FUNCT__ "PetscDTReconstructPoly" 709 /*@ 710 PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 711 712 Not Collective 713 714 Input Arguments: 715 + degree - degree of reconstruction polynomial 716 . nsource - number of source intervals 717 . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 718 . ntarget - number of target intervals 719 - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 720 721 Output Arguments: 722 . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 723 724 Level: advanced 725 726 .seealso: PetscDTLegendreEval() 727 @*/ 728 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 729 { 730 PetscErrorCode ierr; 731 PetscInt i,j,k,*bdegrees,worksize; 732 PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 733 PetscScalar *tau,*work; 734 735 PetscFunctionBegin; 736 PetscValidRealPointer(sourcex,3); 737 PetscValidRealPointer(targetx,5); 738 PetscValidRealPointer(R,6); 739 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); 740 #if defined(PETSC_USE_DEBUG) 741 for (i=0; i<nsource; i++) { 742 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]); 743 } 744 for (i=0; i<ntarget; i++) { 745 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]); 746 } 747 #endif 748 xmin = PetscMin(sourcex[0],targetx[0]); 749 xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 750 center = (xmin + xmax)/2; 751 hscale = (xmax - xmin)/2; 752 worksize = nsource; 753 ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 754 ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 755 for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 756 for (i=0; i<=degree; i++) bdegrees[i] = i+1; 757 ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 758 ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 759 for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 760 ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 761 for (i=0; i<ntarget; i++) { 762 PetscReal rowsum = 0; 763 for (j=0; j<nsource; j++) { 764 PetscReal sum = 0; 765 for (k=0; k<degree+1; k++) { 766 sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 767 } 768 R[i*nsource+j] = sum; 769 rowsum += sum; 770 } 771 for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 772 } 773 ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 774 ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 775 PetscFunctionReturn(0); 776 } 777