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