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