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", (double)quad.points[q*quad.dim+d]);CHKERRQ(ierr); 40 } 41 ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)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 = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr); 151 ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 152 LDZ = N; 153 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 154 PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 155 ierr = PetscFPTrapPop();CHKERRQ(ierr); 156 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 157 158 for (i=0; i<(npoints+1)/2; i++) { 159 PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 160 x[i] = (a+b)/2 - y*(b-a)/2; 161 x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 162 163 w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints])); 164 } 165 ierr = PetscFree2(Z,work);CHKERRQ(ierr); 166 PetscFunctionReturn(0); 167 } 168 169 #undef __FUNCT__ 170 #define __FUNCT__ "PetscDTFactorial_Internal" 171 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 172 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 173 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 174 { 175 PetscReal f = 1.0; 176 PetscInt i; 177 178 PetscFunctionBegin; 179 for (i = 1; i < n+1; ++i) f *= i; 180 *factorial = f; 181 PetscFunctionReturn(0); 182 } 183 184 #undef __FUNCT__ 185 #define __FUNCT__ "PetscDTComputeJacobi" 186 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 187 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 188 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 189 { 190 PetscReal apb, pn1, pn2; 191 PetscInt k; 192 193 PetscFunctionBegin; 194 if (!n) {*P = 1.0; PetscFunctionReturn(0);} 195 if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 196 apb = a + b; 197 pn2 = 1.0; 198 pn1 = 0.5 * (a - b + (apb + 2.0) * x); 199 *P = 0.0; 200 for (k = 2; k < n+1; ++k) { 201 PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 202 PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 203 PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 204 PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 205 206 a2 = a2 / a1; 207 a3 = a3 / a1; 208 a4 = a4 / a1; 209 *P = (a2 + a3 * x) * pn1 - a4 * pn2; 210 pn2 = pn1; 211 pn1 = *P; 212 } 213 PetscFunctionReturn(0); 214 } 215 216 #undef __FUNCT__ 217 #define __FUNCT__ "PetscDTComputeJacobiDerivative" 218 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 219 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 220 { 221 PetscReal nP; 222 PetscErrorCode ierr; 223 224 PetscFunctionBegin; 225 if (!n) {*P = 0.0; PetscFunctionReturn(0);} 226 ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 227 *P = 0.5 * (a + b + n + 1) * nP; 228 PetscFunctionReturn(0); 229 } 230 231 #undef __FUNCT__ 232 #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 233 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 234 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 235 { 236 PetscFunctionBegin; 237 *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 238 *eta = y; 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 244 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 245 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 246 { 247 PetscFunctionBegin; 248 *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 249 *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 250 *zeta = z; 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNCT__ 255 #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 256 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 257 { 258 PetscInt maxIter = 100; 259 PetscReal eps = 1.0e-8; 260 PetscReal a1, a2, a3, a4, a5, a6; 261 PetscInt k; 262 PetscErrorCode ierr; 263 264 PetscFunctionBegin; 265 266 a1 = PetscPowReal(2.0, a+b+1); 267 #if defined(PETSC_HAVE_TGAMMA) 268 a2 = PetscTGamma(a + npoints + 1); 269 a3 = PetscTGamma(b + npoints + 1); 270 a4 = PetscTGamma(a + b + npoints + 1); 271 #else 272 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 273 #endif 274 275 ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 276 a6 = a1 * a2 * a3 / a4 / a5; 277 /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 278 Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 279 for (k = 0; k < npoints; ++k) { 280 PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 281 PetscInt j; 282 283 if (k > 0) r = 0.5 * (r + x[k-1]); 284 for (j = 0; j < maxIter; ++j) { 285 PetscReal s = 0.0, delta, f, fp; 286 PetscInt i; 287 288 for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 289 ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 290 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 291 delta = f / (fp - f * s); 292 r = r - delta; 293 if (PetscAbs(delta) < eps) break; 294 } 295 x[k] = r; 296 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 297 w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 298 } 299 PetscFunctionReturn(0); 300 } 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 304 /*@C 305 PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 306 307 Not Collective 308 309 Input Arguments: 310 + dim - The simplex dimension 311 . order - The quadrature order 312 . a - left end of interval (often-1) 313 - b - right end of interval (often +1) 314 315 Output Arguments: 316 . q - A PetscQuadrature object 317 318 Level: intermediate 319 320 References: 321 Karniadakis and Sherwin. 322 FIAT 323 324 .seealso: PetscDTGaussQuadrature() 325 @*/ 326 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q) 327 { 328 PetscInt npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order; 329 PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 330 PetscInt i, j, k; 331 PetscErrorCode ierr; 332 333 PetscFunctionBegin; 334 if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 335 ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 336 ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 337 switch (dim) { 338 case 0: 339 ierr = PetscFree(x);CHKERRQ(ierr); 340 ierr = PetscFree(w);CHKERRQ(ierr); 341 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 342 ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); 343 x[0] = 0.0; 344 w[0] = 1.0; 345 break; 346 case 1: 347 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr); 348 break; 349 case 2: 350 ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr); 351 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 352 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 353 for (i = 0; i < order; ++i) { 354 for (j = 0; j < order; ++j) { 355 ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr); 356 w[i*order+j] = 0.5 * wx[i] * wy[j]; 357 } 358 } 359 ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 360 break; 361 case 3: 362 ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr); 363 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); 364 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); 365 ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 366 for (i = 0; i < order; ++i) { 367 for (j = 0; j < order; ++j) { 368 for (k = 0; k < order; ++k) { 369 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); 370 w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k]; 371 } 372 } 373 } 374 ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 375 break; 376 default: 377 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 378 } 379 q->dim = dim; 380 q->numPoints = npoints; 381 q->points = x; 382 q->weights = w; 383 PetscFunctionReturn(0); 384 } 385 386 #undef __FUNCT__ 387 #define __FUNCT__ "PetscDTPseudoInverseQR" 388 /* Overwrites A. Can only handle full-rank problems with m>=n 389 * A in column-major format 390 * Ainv in row-major format 391 * tau has length m 392 * worksize must be >= max(1,n) 393 */ 394 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 395 { 396 PetscErrorCode ierr; 397 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 398 PetscScalar *A,*Ainv,*R,*Q,Alpha; 399 400 PetscFunctionBegin; 401 #if defined(PETSC_USE_COMPLEX) 402 { 403 PetscInt i,j; 404 ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 405 for (j=0; j<n; j++) { 406 for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 407 } 408 mstride = m; 409 } 410 #else 411 A = A_in; 412 Ainv = Ainv_out; 413 #endif 414 415 ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 416 ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 417 ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 418 ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 419 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 420 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 421 ierr = PetscFPTrapPop();CHKERRQ(ierr); 422 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 423 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 424 425 /* Extract an explicit representation of Q */ 426 Q = Ainv; 427 ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 428 K = N; /* full rank */ 429 PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 430 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 431 432 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 433 Alpha = 1.0; 434 ldb = lda; 435 PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 436 /* Ainv is Q, overwritten with inverse */ 437 438 #if defined(PETSC_USE_COMPLEX) 439 { 440 PetscInt i; 441 for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 442 ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 443 } 444 #endif 445 PetscFunctionReturn(0); 446 } 447 448 #undef __FUNCT__ 449 #define __FUNCT__ "PetscDTLegendreIntegrate" 450 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 451 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 452 { 453 PetscErrorCode ierr; 454 PetscReal *Bv; 455 PetscInt i,j; 456 457 PetscFunctionBegin; 458 ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 459 /* Point evaluation of L_p on all the source vertices */ 460 ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 461 /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 462 for (i=0; i<ninterval; i++) { 463 for (j=0; j<ndegree; j++) { 464 if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 465 else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 466 } 467 } 468 ierr = PetscFree(Bv);CHKERRQ(ierr); 469 PetscFunctionReturn(0); 470 } 471 472 #undef __FUNCT__ 473 #define __FUNCT__ "PetscDTReconstructPoly" 474 /*@ 475 PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 476 477 Not Collective 478 479 Input Arguments: 480 + degree - degree of reconstruction polynomial 481 . nsource - number of source intervals 482 . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 483 . ntarget - number of target intervals 484 - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 485 486 Output Arguments: 487 . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 488 489 Level: advanced 490 491 .seealso: PetscDTLegendreEval() 492 @*/ 493 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 494 { 495 PetscErrorCode ierr; 496 PetscInt i,j,k,*bdegrees,worksize; 497 PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 498 PetscScalar *tau,*work; 499 500 PetscFunctionBegin; 501 PetscValidRealPointer(sourcex,3); 502 PetscValidRealPointer(targetx,5); 503 PetscValidRealPointer(R,6); 504 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); 505 #if defined(PETSC_USE_DEBUG) 506 for (i=0; i<nsource; i++) { 507 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]); 508 } 509 for (i=0; i<ntarget; i++) { 510 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]); 511 } 512 #endif 513 xmin = PetscMin(sourcex[0],targetx[0]); 514 xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 515 center = (xmin + xmax)/2; 516 hscale = (xmax - xmin)/2; 517 worksize = nsource; 518 ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 519 ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 520 for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 521 for (i=0; i<=degree; i++) bdegrees[i] = i+1; 522 ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 523 ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 524 for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 525 ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 526 for (i=0; i<ntarget; i++) { 527 PetscReal rowsum = 0; 528 for (j=0; j<nsource; j++) { 529 PetscReal sum = 0; 530 for (k=0; k<degree+1; k++) { 531 sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 532 } 533 R[i*nsource+j] = sum; 534 rowsum += sum; 535 } 536 for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 537 } 538 ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 539 ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 540 PetscFunctionReturn(0); 541 } 542