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