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