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 13 #undef __FUNCT__ 14 #define __FUNCT__ "PetscDTLegendreEval" 15 /*@ 16 PetscDTLegendreEval - evaluate Legendre polynomial at points 17 18 Not Collective 19 20 Input Arguments: 21 + npoints - number of spatial points to evaluate at 22 . points - array of locations to evaluate at 23 . ndegree - number of basis degrees to evaluate 24 - degrees - sorted array of degrees to evaluate 25 26 Output Arguments: 27 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 28 . D - row-oriented derivative evaluation matrix (or NULL) 29 - D2 - row-oriented second derivative evaluation matrix (or NULL) 30 31 Level: intermediate 32 33 .seealso: PetscDTGaussQuadrature() 34 @*/ 35 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 36 { 37 PetscInt i,maxdegree; 38 39 PetscFunctionBegin; 40 if (!npoints || !ndegree) PetscFunctionReturn(0); 41 maxdegree = degrees[ndegree-1]; 42 for (i=0; i<npoints; i++) { 43 PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 44 PetscInt j,k; 45 x = points[i]; 46 pm2 = 0; 47 pm1 = 1; 48 pd2 = 0; 49 pd1 = 0; 50 pdd2 = 0; 51 pdd1 = 0; 52 k = 0; 53 if (degrees[k] == 0) { 54 if (B) B[i*ndegree+k] = pm1; 55 if (D) D[i*ndegree+k] = pd1; 56 if (D2) D2[i*ndegree+k] = pdd1; 57 k++; 58 } 59 for (j=1; j<=maxdegree; j++,k++) { 60 PetscReal p,d,dd; 61 p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 62 d = pd2 + (2*j-1)*pm1; 63 dd = pdd2 + (2*j-1)*pd1; 64 pm2 = pm1; 65 pm1 = p; 66 pd2 = pd1; 67 pd1 = d; 68 pdd2 = pdd1; 69 pdd1 = dd; 70 if (degrees[k] == j) { 71 if (B) B[i*ndegree+k] = p; 72 if (D) D[i*ndegree+k] = d; 73 if (D2) D2[i*ndegree+k] = dd; 74 } 75 } 76 } 77 PetscFunctionReturn(0); 78 } 79 80 #undef __FUNCT__ 81 #define __FUNCT__ "PetscDTGaussQuadrature" 82 /*@ 83 PetscDTGaussQuadrature - create Gauss quadrature 84 85 Not Collective 86 87 Input Arguments: 88 + npoints - number of points 89 . a - left end of interval (often-1) 90 - b - right end of interval (often +1) 91 92 Output Arguments: 93 + x - quadrature points 94 - w - quadrature weights 95 96 Level: intermediate 97 98 References: 99 Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 100 101 .seealso: PetscDTLegendreEval() 102 @*/ 103 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 104 { 105 PetscErrorCode ierr; 106 PetscInt i; 107 PetscReal *work; 108 PetscScalar *Z; 109 PetscBLASInt N,LDZ,info; 110 111 PetscFunctionBegin; 112 /* Set up the Golub-Welsch system */ 113 for (i=0; i<npoints; i++) { 114 x[i] = 0; /* diagonal is 0 */ 115 if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 116 } 117 ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 118 ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr); 119 ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 120 LDZ = N; 121 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 122 PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 123 ierr = PetscFPTrapPop();CHKERRQ(ierr); 124 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 125 126 for (i=0; i<(npoints+1)/2; i++) { 127 PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 128 x[i] = (a+b)/2 - y*(b-a)/2; 129 x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 130 131 w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints])); 132 } 133 ierr = PetscFree2(Z,work);CHKERRQ(ierr); 134 PetscFunctionReturn(0); 135 } 136 137 #undef __FUNCT__ 138 #define __FUNCT__ "PetscDTFactorial_Internal" 139 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 140 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 141 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 142 { 143 PetscReal f = 1.0; 144 PetscInt i; 145 146 PetscFunctionBegin; 147 for (i = 1; i < n+1; ++i) f *= i; 148 *factorial = f; 149 PetscFunctionReturn(0); 150 } 151 152 #undef __FUNCT__ 153 #define __FUNCT__ "PetscDTComputeJacobi" 154 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 155 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 156 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 157 { 158 PetscReal apb, pn1, pn2; 159 PetscInt k; 160 161 PetscFunctionBegin; 162 if (!n) {*P = 1.0; PetscFunctionReturn(0);} 163 if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 164 apb = a + b; 165 pn2 = 1.0; 166 pn1 = 0.5 * (a - b + (apb + 2.0) * x); 167 *P = 0.0; 168 for (k = 2; k < n+1; ++k) { 169 PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 170 PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 171 PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 172 PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 173 174 a2 = a2 / a1; 175 a3 = a3 / a1; 176 a4 = a4 / a1; 177 *P = (a2 + a3 * x) * pn1 - a4 * pn2; 178 pn2 = pn1; 179 pn1 = *P; 180 } 181 PetscFunctionReturn(0); 182 } 183 184 #undef __FUNCT__ 185 #define __FUNCT__ "PetscDTComputeJacobiDerivative" 186 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 187 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 188 { 189 PetscReal nP; 190 PetscErrorCode ierr; 191 192 PetscFunctionBegin; 193 if (!n) {*P = 0.0; PetscFunctionReturn(0);} 194 ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 195 *P = 0.5 * (a + b + n + 1) * nP; 196 PetscFunctionReturn(0); 197 } 198 199 #undef __FUNCT__ 200 #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 201 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 202 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 203 { 204 PetscFunctionBegin; 205 *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 206 *eta = y; 207 PetscFunctionReturn(0); 208 } 209 210 #undef __FUNCT__ 211 #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 212 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 213 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 214 { 215 PetscFunctionBegin; 216 *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 217 *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 218 *zeta = z; 219 PetscFunctionReturn(0); 220 } 221 222 #undef __FUNCT__ 223 #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 224 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 225 { 226 PetscInt maxIter = 100; 227 PetscReal eps = 1.0e-8; 228 PetscReal a1, a2, a3, a4, a5, a6; 229 PetscInt k; 230 PetscErrorCode ierr; 231 232 PetscFunctionBegin; 233 234 a1 = pow(2, a+b+1); 235 #if defined(PETSC_HAVE_TGAMMA) 236 a2 = tgamma(a + npoints + 1); 237 a3 = tgamma(b + npoints + 1); 238 a4 = tgamma(a + b + npoints + 1); 239 #else 240 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 241 #endif 242 243 ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 244 a6 = a1 * a2 * a3 / a4 / a5; 245 /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 246 Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 247 for (k = 0; k < npoints; ++k) { 248 PetscReal r = -cos((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 249 PetscInt j; 250 251 if (k > 0) r = 0.5 * (r + x[k-1]); 252 for (j = 0; j < maxIter; ++j) { 253 PetscReal s = 0.0, delta, f, fp; 254 PetscInt i; 255 256 for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 257 ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 258 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 259 delta = f / (fp - f * s); 260 r = r - delta; 261 if (fabs(delta) < eps) break; 262 } 263 x[k] = r; 264 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 265 w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 266 } 267 PetscFunctionReturn(0); 268 } 269 270 #undef __FUNCT__ 271 #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 272 /*@C 273 PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 274 275 Not Collective 276 277 Input Arguments: 278 + dim - The simplex dimension 279 . npoints - number of points 280 . a - left end of interval (often-1) 281 - b - right end of interval (often +1) 282 283 Output Arguments: 284 + points - quadrature points 285 - weights - quadrature weights 286 287 Level: intermediate 288 289 References: 290 Karniadakis and Sherwin. 291 FIAT 292 293 .seealso: PetscDTGaussQuadrature() 294 @*/ 295 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscReal *points[], PetscReal *weights[]) 296 { 297 PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 298 PetscInt i, j, k; 299 PetscErrorCode ierr; 300 301 PetscFunctionBegin; 302 if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 303 switch (dim) { 304 case 1: 305 ierr = PetscMalloc(npoints * sizeof(PetscReal), &x);CHKERRQ(ierr); 306 ierr = PetscMalloc(npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 307 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, w);CHKERRQ(ierr); 308 break; 309 case 2: 310 ierr = PetscMalloc(npoints*npoints*2 * sizeof(PetscReal), &x);CHKERRQ(ierr); 311 ierr = PetscMalloc(npoints*npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 312 ierr = PetscMalloc4(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy);CHKERRQ(ierr); 313 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 314 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 315 for (i = 0; i < npoints; ++i) { 316 for (j = 0; j < npoints; ++j) { 317 ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 318 w[i*npoints+j] = 0.5 * wx[i] * wy[j]; 319 } 320 } 321 ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 322 break; 323 case 3: 324 ierr = PetscMalloc(npoints*npoints*3 * sizeof(PetscReal), &x);CHKERRQ(ierr); 325 ierr = PetscMalloc(npoints*npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 326 ierr = PetscMalloc6(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy,npoints,PetscReal,&pz,npoints,PetscReal,&wz);CHKERRQ(ierr); 327 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 328 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 329 ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 330 for (i = 0; i < npoints; ++i) { 331 for (j = 0; j < npoints; ++j) { 332 for (k = 0; k < npoints; ++k) { 333 ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*npoints+j)*npoints+k)*3+0], &x[((i*npoints+j)*npoints+k)*3+1], &x[((i*npoints+j)*npoints+k)*3+2]);CHKERRQ(ierr); 334 w[(i*npoints+j)*npoints+k] = 0.125 * wx[i] * wy[j] * wz[k]; 335 } 336 } 337 } 338 ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 339 break; 340 default: 341 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 342 } 343 if (points) *points = x; 344 if (weights) *weights = w; 345 PetscFunctionReturn(0); 346 } 347 348 #undef __FUNCT__ 349 #define __FUNCT__ "PetscDTPseudoInverseQR" 350 /* Overwrites A. Can only handle full-rank problems with m>=n 351 * A in column-major format 352 * Ainv in row-major format 353 * tau has length m 354 * worksize must be >= max(1,n) 355 */ 356 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 357 { 358 PetscErrorCode ierr; 359 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 360 PetscScalar *A,*Ainv,*R,*Q,Alpha; 361 362 PetscFunctionBegin; 363 #if defined(PETSC_USE_COMPLEX) 364 { 365 PetscInt i,j; 366 ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr); 367 for (j=0; j<n; j++) { 368 for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 369 } 370 mstride = m; 371 } 372 #else 373 A = A_in; 374 Ainv = Ainv_out; 375 #endif 376 377 ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 378 ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 379 ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 380 ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 381 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 382 LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info); 383 ierr = PetscFPTrapPop();CHKERRQ(ierr); 384 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 385 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 386 387 /* Extract an explicit representation of Q */ 388 Q = Ainv; 389 ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 390 K = N; /* full rank */ 391 LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info); 392 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 393 394 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 395 Alpha = 1.0; 396 ldb = lda; 397 BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb); 398 /* Ainv is Q, overwritten with inverse */ 399 400 #if defined(PETSC_USE_COMPLEX) 401 { 402 PetscInt i; 403 for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 404 ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 405 } 406 #endif 407 PetscFunctionReturn(0); 408 } 409 410 #undef __FUNCT__ 411 #define __FUNCT__ "PetscDTLegendreIntegrate" 412 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 413 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 414 { 415 PetscErrorCode ierr; 416 PetscReal *Bv; 417 PetscInt i,j; 418 419 PetscFunctionBegin; 420 ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr); 421 /* Point evaluation of L_p on all the source vertices */ 422 ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 423 /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 424 for (i=0; i<ninterval; i++) { 425 for (j=0; j<ndegree; j++) { 426 if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 427 else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 428 } 429 } 430 ierr = PetscFree(Bv);CHKERRQ(ierr); 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "PetscDTReconstructPoly" 436 /*@ 437 PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 438 439 Not Collective 440 441 Input Arguments: 442 + degree - degree of reconstruction polynomial 443 . nsource - number of source intervals 444 . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 445 . ntarget - number of target intervals 446 - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 447 448 Output Arguments: 449 . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 450 451 Level: advanced 452 453 .seealso: PetscDTLegendreEval() 454 @*/ 455 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 456 { 457 PetscErrorCode ierr; 458 PetscInt i,j,k,*bdegrees,worksize; 459 PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 460 PetscScalar *tau,*work; 461 462 PetscFunctionBegin; 463 PetscValidRealPointer(sourcex,3); 464 PetscValidRealPointer(targetx,5); 465 PetscValidRealPointer(R,6); 466 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); 467 #if defined(PETSC_USE_DEBUG) 468 for (i=0; i<nsource; i++) { 469 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]); 470 } 471 for (i=0; i<ntarget; i++) { 472 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]); 473 } 474 #endif 475 xmin = PetscMin(sourcex[0],targetx[0]); 476 xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 477 center = (xmin + xmax)/2; 478 hscale = (xmax - xmin)/2; 479 worksize = nsource; 480 ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr); 481 ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr); 482 for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 483 for (i=0; i<=degree; i++) bdegrees[i] = i+1; 484 ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 485 ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 486 for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 487 ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 488 for (i=0; i<ntarget; i++) { 489 PetscReal rowsum = 0; 490 for (j=0; j<nsource; j++) { 491 PetscReal sum = 0; 492 for (k=0; k<degree+1; k++) { 493 sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 494 } 495 R[i*nsource+j] = sum; 496 rowsum += sum; 497 } 498 for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 499 } 500 ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 501 ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 502 PetscFunctionReturn(0); 503 } 504 505 /* Basis Jet Tabulation 506 507 We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We 508 follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can 509 be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis 510 as a prime basis. 511 512 \psi_i = \sum_k \alpha_{ki} \phi_k 513 514 Our nodal basis is defined in terms of the dual basis $n_j$ 515 516 n_j \cdot \psi_i = \delta_{ji} 517 518 and we may act on the first equation to obtain 519 520 n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k 521 \delta_{ji} = \sum_k \alpha_{ki} V_{jk} 522 I = V \alpha 523 524 so the coefficients of the nodal basis in the prime basis are 525 526 \alpha = V^{-1} 527 528 We will define the dual basis vectors $n_j$ using a quadrature rule. 529 530 Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials 531 (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can 532 be implemented exactly as in FIAT using functionals $L_j$. 533 534 I will have to count the degrees correctly for the Legendre product when we are on simplices. 535 536 We will have three objects: 537 - Space, P: this just need point evaluation I think 538 - Dual Space, P'+K: This looks like a set of functionals that can act on members of P, each n is defined by a Q 539 - FEM: This keeps {P, P', Q} 540 */ 541 #include <petsc-private/petscfeimpl.h> 542 543 PetscInt PETSCSPACE_CLASSID = 0; 544 545 PetscFunctionList PetscSpaceList = NULL; 546 PetscBool PetscSpaceRegisterAllCalled = PETSC_FALSE; 547 548 #undef __FUNCT__ 549 #define __FUNCT__ "PetscSpaceRegister" 550 /*@C 551 PetscSpaceRegister - Adds a new PetscSpace implementation 552 553 Not Collective 554 555 Input Parameters: 556 + name - The name of a new user-defined creation routine 557 - create_func - The creation routine itself 558 559 Notes: 560 PetscSpaceRegister() may be called multiple times to add several user-defined PetscSpaces 561 562 Sample usage: 563 .vb 564 PetscSpaceRegister("my_space", MyPetscSpaceCreate); 565 .ve 566 567 Then, your PetscSpace type can be chosen with the procedural interface via 568 .vb 569 PetscSpaceCreate(MPI_Comm, PetscSpace *); 570 PetscSpaceSetType(PetscSpace, "my_space"); 571 .ve 572 or at runtime via the option 573 .vb 574 -petscspace_type my_space 575 .ve 576 577 Level: advanced 578 579 .keywords: PetscSpace, register 580 .seealso: PetscSpaceRegisterAll(), PetscSpaceRegisterDestroy() 581 582 @*/ 583 PetscErrorCode PetscSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscSpace)) 584 { 585 PetscErrorCode ierr; 586 587 PetscFunctionBegin; 588 ierr = PetscFunctionListAdd(&PetscSpaceList, sname, function);CHKERRQ(ierr); 589 PetscFunctionReturn(0); 590 } 591 592 #undef __FUNCT__ 593 #define __FUNCT__ "PetscSpaceSetType" 594 /*@C 595 PetscSpaceSetType - Builds a particular PetscSpace 596 597 Collective on PetscSpace 598 599 Input Parameters: 600 + sp - The PetscSpace object 601 - name - The kind of space 602 603 Options Database Key: 604 . -petscspace_type <type> - Sets the PetscSpace type; use -help for a list of available types 605 606 Level: intermediate 607 608 .keywords: PetscSpace, set, type 609 .seealso: PetscSpaceGetType(), PetscSpaceCreate() 610 @*/ 611 PetscErrorCode PetscSpaceSetType(PetscSpace sp, PetscSpaceType name) 612 { 613 PetscErrorCode (*r)(PetscSpace); 614 PetscBool match; 615 PetscErrorCode ierr; 616 617 PetscFunctionBegin; 618 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 619 ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 620 if (match) PetscFunctionReturn(0); 621 622 if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 623 ierr = PetscFunctionListFind(PetscSpaceList, name, &r);CHKERRQ(ierr); 624 if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSpace type: %s", name); 625 626 if (sp->ops->destroy) { 627 ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 628 sp->ops->destroy = NULL; 629 } 630 ierr = (*r)(sp);CHKERRQ(ierr); 631 ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 632 PetscFunctionReturn(0); 633 } 634 635 #undef __FUNCT__ 636 #define __FUNCT__ "PetscSpaceGetType" 637 /*@C 638 PetscSpaceGetType - Gets the PetscSpace type name (as a string) from the object. 639 640 Not Collective 641 642 Input Parameter: 643 . dm - The PetscSpace 644 645 Output Parameter: 646 . name - The PetscSpace type name 647 648 Level: intermediate 649 650 .keywords: PetscSpace, get, type, name 651 .seealso: PetscSpaceSetType(), PetscSpaceCreate() 652 @*/ 653 PetscErrorCode PetscSpaceGetType(PetscSpace sp, PetscSpaceType *name) 654 { 655 PetscErrorCode ierr; 656 657 PetscFunctionBegin; 658 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 659 PetscValidCharPointer(name, 2); 660 if (!PetscSpaceRegisterAllCalled) { 661 ierr = PetscSpaceRegisterAll();CHKERRQ(ierr); 662 } 663 *name = ((PetscObject) sp)->type_name; 664 PetscFunctionReturn(0); 665 } 666 667 #undef __FUNCT__ 668 #define __FUNCT__ "PetscSpaceView" 669 /*@C 670 PetscSpaceView - Views a PetscSpace 671 672 Collective on PetscSpace 673 674 Input Parameter: 675 + sp - the PetscSpace object to view 676 - v - the viewer 677 678 Level: developer 679 680 .seealso PetscSpaceDestroy() 681 @*/ 682 PetscErrorCode PetscSpaceView(PetscSpace sp, PetscViewer v) 683 { 684 PetscErrorCode ierr; 685 686 PetscFunctionBegin; 687 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 688 if (!v) { 689 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr); 690 } 691 if (sp->ops->view) { 692 ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr); 693 } 694 PetscFunctionReturn(0); 695 } 696 697 #undef __FUNCT__ 698 #define __FUNCT__ "PetscSpaceViewFromOptions" 699 /* 700 PetscSpaceViewFromOptions - Processes command line options to determine if/how a PetscSpace is to be viewed. 701 702 Collective on PetscSpace 703 704 Input Parameters: 705 + sp - the PetscSpace 706 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd' 707 - optionname - option to activate viewing 708 709 Level: intermediate 710 711 .keywords: PetscSpace, view, options, database 712 .seealso: VecViewFromOptions(), MatViewFromOptions() 713 */ 714 PetscErrorCode PetscSpaceViewFromOptions(PetscSpace sp, const char prefix[], const char optionname[]) 715 { 716 PetscViewer viewer; 717 PetscViewerFormat format; 718 PetscBool flg; 719 PetscErrorCode ierr; 720 721 PetscFunctionBegin; 722 if (prefix) { 723 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 724 } else { 725 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 726 } 727 if (flg) { 728 ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 729 ierr = PetscSpaceView(sp, viewer);CHKERRQ(ierr); 730 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 731 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 732 } 733 PetscFunctionReturn(0); 734 } 735 736 #undef __FUNCT__ 737 #define __FUNCT__ "PetscSpaceSetFromOptions" 738 /*@ 739 PetscSpaceSetFromOptions - sets parameters in a PetscSpace from the options database 740 741 Collective on PetscSpace 742 743 Input Parameter: 744 . sp - the PetscSpace object to set options for 745 746 Options Database: 747 . -petscspace_order the approximation order of the space 748 749 Level: developer 750 751 .seealso PetscSpaceView() 752 @*/ 753 PetscErrorCode PetscSpaceSetFromOptions(PetscSpace sp) 754 { 755 const char *defaultType; 756 char typename[256]; 757 PetscBool flg; 758 PetscErrorCode ierr; 759 760 PetscFunctionBegin; 761 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 762 if (!((PetscObject) sp)->type_name) { 763 defaultType = PETSCSPACEPOLYNOMIAL; 764 } else { 765 defaultType = ((PetscObject) sp)->type_name; 766 } 767 if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 768 769 ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 770 ierr = PetscOptionsList("-petscspace_type", "Linear space", "PetscSpaceSetType", PetscSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr); 771 if (flg) { 772 ierr = PetscSpaceSetType(sp, typename);CHKERRQ(ierr); 773 } else if (!((PetscObject) sp)->type_name) { 774 ierr = PetscSpaceSetType(sp, defaultType);CHKERRQ(ierr); 775 } 776 ierr = PetscOptionsInt("-petscspace_order", "The approximation order", "PetscSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr); 777 if (sp->ops->setfromoptions) { 778 ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr); 779 } 780 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 781 ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr); 782 ierr = PetscOptionsEnd();CHKERRQ(ierr); 783 ierr = PetscSpaceViewFromOptions(sp, NULL, "-petscspace_view");CHKERRQ(ierr); 784 PetscFunctionReturn(0); 785 } 786 787 #undef __FUNCT__ 788 #define __FUNCT__ "PetscSpaceSetUp" 789 /*@C 790 PetscSpaceSetUp - Construct data structures for the PetscSpace 791 792 Collective on PetscSpace 793 794 Input Parameter: 795 . sp - the PetscSpace object to setup 796 797 Level: developer 798 799 .seealso PetscSpaceView(), PetscSpaceDestroy() 800 @*/ 801 PetscErrorCode PetscSpaceSetUp(PetscSpace sp) 802 { 803 PetscErrorCode ierr; 804 805 PetscFunctionBegin; 806 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 807 if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 808 PetscFunctionReturn(0); 809 } 810 811 #undef __FUNCT__ 812 #define __FUNCT__ "PetscSpaceDestroy" 813 /*@ 814 PetscSpaceDestroy - Destroys a PetscSpace object 815 816 Collective on PetscSpace 817 818 Input Parameter: 819 . sp - the PetscSpace object to destroy 820 821 Level: developer 822 823 .seealso PetscSpaceView() 824 @*/ 825 PetscErrorCode PetscSpaceDestroy(PetscSpace *sp) 826 { 827 PetscErrorCode ierr; 828 829 PetscFunctionBegin; 830 if (!*sp) PetscFunctionReturn(0); 831 PetscValidHeaderSpecific((*sp), PETSCSPACE_CLASSID, 1); 832 833 if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 834 ((PetscObject) (*sp))->refct = 0; 835 /* if memory was published with AMS then destroy it */ 836 ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr); 837 838 ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr); 839 840 ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr); 841 ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 842 PetscFunctionReturn(0); 843 } 844 845 #undef __FUNCT__ 846 #define __FUNCT__ "PetscSpaceCreate" 847 /*@ 848 PetscSpaceCreate - Creates an empty PetscSpace object. The type can then be set with PetscSpaceSetType(). 849 850 Collective on MPI_Comm 851 852 Input Parameter: 853 . comm - The communicator for the PetscSpace object 854 855 Output Parameter: 856 . sp - The PetscSpace object 857 858 Level: beginner 859 860 .seealso: PetscSpaceSetType(), PETSCSPACEPOLYNOMIAL 861 @*/ 862 PetscErrorCode PetscSpaceCreate(MPI_Comm comm, PetscSpace *sp) 863 { 864 PetscSpace s; 865 PetscErrorCode ierr; 866 867 PetscFunctionBegin; 868 PetscValidPointer(sp, 2); 869 *sp = NULL; 870 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 871 ierr = PetscFEInitializePackage();CHKERRQ(ierr); 872 #endif 873 874 ierr = PetscHeaderCreate(s, _p_PetscSpace, struct _PetscSpaceOps, PETSCSPACE_CLASSID, "PetscSpace", "Linear Space", "PetscSpace", comm, PetscSpaceDestroy, PetscSpaceView);CHKERRQ(ierr); 875 ierr = PetscMemzero(s->ops, sizeof(struct _PetscSpaceOps));CHKERRQ(ierr); 876 877 s->order = 0; 878 ierr = DMShellCreate(comm, &s->dm);CHKERRQ(ierr); 879 880 *sp = s; 881 PetscFunctionReturn(0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "PetscSpaceGetDimension" 886 /* Dimension of the space, i.e. number of basis vectors */ 887 PetscErrorCode PetscSpaceGetDimension(PetscSpace sp, PetscInt *dim) 888 { 889 PetscErrorCode ierr; 890 891 PetscFunctionBegin; 892 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 893 PetscValidPointer(dim, 2); 894 *dim = 0; 895 if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 896 PetscFunctionReturn(0); 897 } 898 899 #undef __FUNCT__ 900 #define __FUNCT__ "PetscSpaceGetOrder" 901 PetscErrorCode PetscSpaceGetOrder(PetscSpace sp, PetscInt *order) 902 { 903 PetscFunctionBegin; 904 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 905 PetscValidPointer(order, 2); 906 *order = sp->order; 907 PetscFunctionReturn(0); 908 } 909 910 #undef __FUNCT__ 911 #define __FUNCT__ "PetscSpaceSetOrder" 912 PetscErrorCode PetscSpaceSetOrder(PetscSpace sp, PetscInt order) 913 { 914 PetscFunctionBegin; 915 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 916 sp->order = order; 917 PetscFunctionReturn(0); 918 } 919 920 #undef __FUNCT__ 921 #define __FUNCT__ "PetscSpaceEvaluate" 922 PetscErrorCode PetscSpaceEvaluate(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 923 { 924 PetscErrorCode ierr; 925 926 PetscFunctionBegin; 927 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 928 PetscValidPointer(points, 3); 929 if (B) PetscValidPointer(B, 4); 930 if (D) PetscValidPointer(D, 5); 931 if (H) PetscValidPointer(H, 6); 932 if (sp->ops->evaluate) {ierr = (*sp->ops->evaluate)(sp, npoints, points, B, D, H);CHKERRQ(ierr);} 933 PetscFunctionReturn(0); 934 } 935 936 #undef __FUNCT__ 937 #define __FUNCT__ "PetscSpaceSetFromOptions_Polynomial" 938 PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscSpace sp) 939 { 940 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 941 PetscErrorCode ierr; 942 943 PetscFunctionBegin; 944 ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 945 ierr = PetscOptionsInt("-petscspace_num_variables", "The number of different variables, e.g. x and y", "PetscSpacePolynomialSetNumVariables", poly->numVariables, &poly->numVariables, NULL);CHKERRQ(ierr); 946 ierr = PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);CHKERRQ(ierr); 947 ierr = PetscOptionsEnd();CHKERRQ(ierr); 948 PetscFunctionReturn(0); 949 } 950 951 #undef __FUNCT__ 952 #define __FUNCT__ "PetscSpacePolynomialView_Ascii" 953 PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer viewer) 954 { 955 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 956 PetscViewerFormat format; 957 PetscErrorCode ierr; 958 959 PetscFunctionBegin; 960 ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 961 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 962 ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr); 963 } else { 964 ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr); 965 } 966 PetscFunctionReturn(0); 967 } 968 969 #undef __FUNCT__ 970 #define __FUNCT__ "PetscSpaceView_Polynomial" 971 PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer) 972 { 973 PetscBool iascii; 974 PetscErrorCode ierr; 975 976 PetscFunctionBegin; 977 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 978 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 979 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 980 if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);} 981 PetscFunctionReturn(0); 982 } 983 984 #undef __FUNCT__ 985 #define __FUNCT__ "PetscSpaceSetUp_Polynomial" 986 PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp) 987 { 988 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 989 PetscInt ndegree = sp->order+1; 990 PetscInt dim = poly->numVariables; 991 PetscInt pdim, deg; 992 PetscErrorCode ierr; 993 994 PetscFunctionBegin; 995 ierr = PetscMalloc(ndegree * sizeof(PetscInt), &poly->degrees);CHKERRQ(ierr); 996 for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg; 997 ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 998 /* Create A (pdim x ndegree * dim) */ 999 ierr = PetscMalloc(pdim*ndegree*dim * sizeof(PetscReal), &poly->A);CHKERRQ(ierr); 1000 ierr = PetscMemzero(poly->A, pdim*ndegree*dim * sizeof(PetscReal));CHKERRQ(ierr); 1001 /* Hardcode P_1: Here we need a way to iterate through the basis */ 1002 poly->A[(0*ndegree + 0)*dim + 0] = 1.0; /* 1 */ 1003 poly->A[(1*ndegree + 1)*dim + 0] = 1.0; /* x */ 1004 poly->A[(2*ndegree + 1)*dim + 1] = 1.0; /* y */ 1005 PetscFunctionReturn(0); 1006 } 1007 1008 #undef __FUNCT__ 1009 #define __FUNCT__ "PetscSpaceDestroy_Polynomial" 1010 PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp) 1011 { 1012 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1013 PetscErrorCode ierr; 1014 1015 PetscFunctionBegin; 1016 ierr = PetscFree(poly->A);CHKERRQ(ierr); 1017 ierr = PetscFree(poly->degrees);CHKERRQ(ierr); 1018 ierr = PetscFree(poly);CHKERRQ(ierr); 1019 PetscFunctionReturn(0); 1020 } 1021 1022 #undef __FUNCT__ 1023 #define __FUNCT__ "PetscSpaceGetDimension_Polynomial" 1024 PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim) 1025 { 1026 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1027 PetscInt deg = sp->order; 1028 PetscInt n = poly->numVariables, i; 1029 PetscReal D = 1.0; 1030 1031 PetscFunctionBegin; 1032 for (i = 1; i <= n; ++i) { 1033 D *= ((PetscReal) (deg+i))/i; 1034 } 1035 *dim = (PetscInt) (D + 0.5); 1036 PetscFunctionReturn(0); 1037 } 1038 1039 #undef __FUNCT__ 1040 #define __FUNCT__ "PetscSpaceEvaluate_Polynomial" 1041 PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 1042 { 1043 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1044 DM dm = sp->dm; 1045 PetscInt ndegree = sp->order+1; 1046 PetscInt *degrees = poly->degrees; 1047 PetscInt dim = poly->numVariables; 1048 PetscReal *A = poly->A; 1049 PetscReal *lpoints, *tmp, *LB, *LD, *LH; 1050 PetscInt pdim, d, i, p, deg; 1051 PetscErrorCode ierr; 1052 1053 PetscFunctionBegin; 1054 ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 1055 ierr = DMGetWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr); 1056 ierr = DMGetWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr); 1057 if (B) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);} 1058 if (D) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);} 1059 if (H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);} 1060 for (d = 0; d < dim; ++d) { 1061 for (p = 0; p < npoints; ++p) { 1062 lpoints[p] = points[p*dim+d]; 1063 } 1064 ierr = PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);CHKERRQ(ierr); 1065 /* LB, LD, LH (ndegree * dim x npoints) */ 1066 for (deg = 0; deg < ndegree; ++deg) { 1067 for (p = 0; p < npoints; ++p) { 1068 if (B) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg]; 1069 if (D) LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg]; 1070 if (H) LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg]; 1071 } 1072 } 1073 } 1074 /* Multiply by A (pdim x ndegree * dim) */ 1075 if (B) { 1076 for (p = 0; p < npoints; ++p) { 1077 for (i = 0; i < pdim; ++i) { 1078 B[p*pdim + i] = 0.0; 1079 for (deg = 0; deg < ndegree; ++deg) { 1080 for (d = 0; d < dim; ++d) { 1081 B[p*pdim + i] += A[(i*ndegree + deg)*dim + d] * LB[(deg*dim + d)*npoints + p]; 1082 } 1083 } 1084 } 1085 } 1086 } 1087 if (D) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code first derivatives"); 1088 if (H) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code second derivatives"); 1089 if (B) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);} 1090 if (D) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);} 1091 if (H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);} 1092 ierr = DMRestoreWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr); 1093 ierr = DMRestoreWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr); 1094 PetscFunctionReturn(0); 1095 } 1096 1097 #undef __FUNCT__ 1098 #define __FUNCT__ "PetscSpaceInitialize_Polynomial" 1099 PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp) 1100 { 1101 PetscFunctionBegin; 1102 sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial; 1103 sp->ops->setup = PetscSpaceSetUp_Polynomial; 1104 sp->ops->view = PetscSpaceView_Polynomial; 1105 sp->ops->destroy = PetscSpaceDestroy_Polynomial; 1106 sp->ops->getdimension = PetscSpaceGetDimension_Polynomial; 1107 sp->ops->evaluate = PetscSpaceEvaluate_Polynomial; 1108 PetscFunctionReturn(0); 1109 } 1110 1111 /*MC 1112 PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of linear polynomials. 1113 1114 Level: intermediate 1115 1116 .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType() 1117 M*/ 1118 1119 #undef __FUNCT__ 1120 #define __FUNCT__ "PetscSpaceCreate_Polynomial" 1121 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp) 1122 { 1123 PetscSpace_Poly *poly; 1124 PetscErrorCode ierr; 1125 1126 PetscFunctionBegin; 1127 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1128 ierr = PetscNewLog(sp, PetscSpace_Poly, &poly);CHKERRQ(ierr); 1129 sp->data = poly; 1130 1131 poly->numVariables = 0; 1132 poly->symmetric = PETSC_FALSE; 1133 poly->degrees = NULL; 1134 poly->A = NULL; 1135 1136 ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr); 1137 PetscFunctionReturn(0); 1138 } 1139 1140 #undef __FUNCT__ 1141 #define __FUNCT__ "PetscSpacePolynomialSetSymmetric" 1142 PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool sym) 1143 { 1144 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1145 1146 PetscFunctionBegin; 1147 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1148 poly->symmetric = sym; 1149 PetscFunctionReturn(0); 1150 } 1151 1152 #undef __FUNCT__ 1153 #define __FUNCT__ "PetscSpacePolynomialGetSymmetric" 1154 PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace sp, PetscBool *sym) 1155 { 1156 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1157 1158 PetscFunctionBegin; 1159 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1160 PetscValidPointer(sym, 2); 1161 *sym = poly->symmetric; 1162 PetscFunctionReturn(0); 1163 } 1164 1165 #undef __FUNCT__ 1166 #define __FUNCT__ "PetscSpacePolynomialSetNumVariables" 1167 PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace sp, PetscInt n) 1168 { 1169 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1170 1171 PetscFunctionBegin; 1172 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1173 poly->numVariables = n; 1174 PetscFunctionReturn(0); 1175 } 1176 1177 #undef __FUNCT__ 1178 #define __FUNCT__ "PetscSpacePolynomialGetNumVariables" 1179 PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace sp, PetscInt *n) 1180 { 1181 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1182 1183 PetscFunctionBegin; 1184 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1185 PetscValidPointer(n, 2); 1186 *n = poly->numVariables; 1187 PetscFunctionReturn(0); 1188 } 1189 1190 1191 PetscInt PETSCDUALSPACE_CLASSID = 0; 1192 1193 PetscFunctionList PetscDualSpaceList = NULL; 1194 PetscBool PetscDualSpaceRegisterAllCalled = PETSC_FALSE; 1195 1196 #undef __FUNCT__ 1197 #define __FUNCT__ "PetscDualSpaceRegister" 1198 /*@C 1199 PetscDualSpaceRegister - Adds a new PetscDualSpace implementation 1200 1201 Not Collective 1202 1203 Input Parameters: 1204 + name - The name of a new user-defined creation routine 1205 - create_func - The creation routine itself 1206 1207 Notes: 1208 PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces 1209 1210 Sample usage: 1211 .vb 1212 PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate); 1213 .ve 1214 1215 Then, your PetscDualSpace type can be chosen with the procedural interface via 1216 .vb 1217 PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 1218 PetscDualSpaceSetType(PetscDualSpace, "my_dual_space"); 1219 .ve 1220 or at runtime via the option 1221 .vb 1222 -petscdualspace_type my_dual_space 1223 .ve 1224 1225 Level: advanced 1226 1227 .keywords: PetscDualSpace, register 1228 .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy() 1229 1230 @*/ 1231 PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace)) 1232 { 1233 PetscErrorCode ierr; 1234 1235 PetscFunctionBegin; 1236 ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr); 1237 PetscFunctionReturn(0); 1238 } 1239 1240 #undef __FUNCT__ 1241 #define __FUNCT__ "PetscDualSpaceSetType" 1242 /*@C 1243 PetscDualSpaceSetType - Builds a particular PetscDualSpace 1244 1245 Collective on PetscDualSpace 1246 1247 Input Parameters: 1248 + sp - The PetscDualSpace object 1249 - name - The kind of space 1250 1251 Options Database Key: 1252 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types 1253 1254 Level: intermediate 1255 1256 .keywords: PetscDualSpace, set, type 1257 .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate() 1258 @*/ 1259 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name) 1260 { 1261 PetscErrorCode (*r)(PetscDualSpace); 1262 PetscBool match; 1263 PetscErrorCode ierr; 1264 1265 PetscFunctionBegin; 1266 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1267 ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 1268 if (match) PetscFunctionReturn(0); 1269 1270 if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);} 1271 ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr); 1272 if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name); 1273 1274 if (sp->ops->destroy) { 1275 ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 1276 sp->ops->destroy = NULL; 1277 } 1278 ierr = (*r)(sp);CHKERRQ(ierr); 1279 ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 1280 PetscFunctionReturn(0); 1281 } 1282 1283 #undef __FUNCT__ 1284 #define __FUNCT__ "PetscDualSpaceGetType" 1285 /*@C 1286 PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object. 1287 1288 Not Collective 1289 1290 Input Parameter: 1291 . dm - The PetscDualSpace 1292 1293 Output Parameter: 1294 . name - The PetscDualSpace type name 1295 1296 Level: intermediate 1297 1298 .keywords: PetscDualSpace, get, type, name 1299 .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate() 1300 @*/ 1301 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name) 1302 { 1303 PetscErrorCode ierr; 1304 1305 PetscFunctionBegin; 1306 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1307 PetscValidCharPointer(name, 2); 1308 if (!PetscDualSpaceRegisterAllCalled) { 1309 ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr); 1310 } 1311 *name = ((PetscObject) sp)->type_name; 1312 PetscFunctionReturn(0); 1313 } 1314 1315 #undef __FUNCT__ 1316 #define __FUNCT__ "PetscDualSpaceView" 1317 /*@C 1318 PetscDualSpaceView - Views a PetscDualSpace 1319 1320 Collective on PetscDualSpace 1321 1322 Input Parameter: 1323 + sp - the PetscDualSpace object to view 1324 - v - the viewer 1325 1326 Level: developer 1327 1328 .seealso PetscDualSpaceDestroy() 1329 @*/ 1330 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v) 1331 { 1332 PetscErrorCode ierr; 1333 1334 PetscFunctionBegin; 1335 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1336 if (!v) { 1337 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr); 1338 } 1339 if (sp->ops->view) { 1340 ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr); 1341 } 1342 PetscFunctionReturn(0); 1343 } 1344 1345 #undef __FUNCT__ 1346 #define __FUNCT__ "PetscDualSpaceViewFromOptions" 1347 /* 1348 PetscDualSpaceViewFromOptions - Processes command line options to determine if/how a PetscDualSpace is to be viewed. 1349 1350 Collective on PetscDualSpace 1351 1352 Input Parameters: 1353 + sp - the PetscDualSpace 1354 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd' 1355 - optionname - option to activate viewing 1356 1357 Level: intermediate 1358 1359 .keywords: PetscDualSpace, view, options, database 1360 .seealso: VecViewFromOptions(), MatViewFromOptions() 1361 */ 1362 PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace sp, const char prefix[], const char optionname[]) 1363 { 1364 PetscViewer viewer; 1365 PetscViewerFormat format; 1366 PetscBool flg; 1367 PetscErrorCode ierr; 1368 1369 PetscFunctionBegin; 1370 if (prefix) { 1371 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 1372 } else { 1373 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 1374 } 1375 if (flg) { 1376 ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 1377 ierr = PetscDualSpaceView(sp, viewer);CHKERRQ(ierr); 1378 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1379 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1380 } 1381 PetscFunctionReturn(0); 1382 } 1383 1384 #undef __FUNCT__ 1385 #define __FUNCT__ "PetscDualSpaceSetFromOptions" 1386 /*@ 1387 PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database 1388 1389 Collective on PetscDualSpace 1390 1391 Input Parameter: 1392 . sp - the PetscDualSpace object to set options for 1393 1394 Options Database: 1395 . -petscspace_order the approximation order of the space 1396 1397 Level: developer 1398 1399 .seealso PetscDualSpaceView() 1400 @*/ 1401 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp) 1402 { 1403 const char *defaultType; 1404 char typename[256]; 1405 PetscBool flg; 1406 PetscErrorCode ierr; 1407 1408 PetscFunctionBegin; 1409 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1410 if (!((PetscObject) sp)->type_name) { 1411 defaultType = PETSCDUALSPACELAGRANGE; 1412 } else { 1413 defaultType = ((PetscObject) sp)->type_name; 1414 } 1415 if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 1416 1417 ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 1418 ierr = PetscOptionsList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr); 1419 if (flg) { 1420 ierr = PetscDualSpaceSetType(sp, typename);CHKERRQ(ierr); 1421 } else if (!((PetscObject) sp)->type_name) { 1422 ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr); 1423 } 1424 ierr = PetscOptionsInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr); 1425 if (sp->ops->setfromoptions) { 1426 ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr); 1427 } 1428 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 1429 ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr); 1430 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1431 ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr); 1432 PetscFunctionReturn(0); 1433 } 1434 1435 #undef __FUNCT__ 1436 #define __FUNCT__ "PetscDualSpaceSetUp" 1437 /*@C 1438 PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace 1439 1440 Collective on PetscDualSpace 1441 1442 Input Parameter: 1443 . sp - the PetscDualSpace object to setup 1444 1445 Level: developer 1446 1447 .seealso PetscDualSpaceView(), PetscDualSpaceDestroy() 1448 @*/ 1449 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp) 1450 { 1451 PetscErrorCode ierr; 1452 1453 PetscFunctionBegin; 1454 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1455 if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 1456 PetscFunctionReturn(0); 1457 } 1458 1459 #undef __FUNCT__ 1460 #define __FUNCT__ "PetscDualSpaceDestroy" 1461 /*@ 1462 PetscDualSpaceDestroy - Destroys a PetscDualSpace object 1463 1464 Collective on PetscDualSpace 1465 1466 Input Parameter: 1467 . sp - the PetscDualSpace object to destroy 1468 1469 Level: developer 1470 1471 .seealso PetscDualSpaceView() 1472 @*/ 1473 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp) 1474 { 1475 PetscInt dim, f; 1476 PetscErrorCode ierr; 1477 1478 PetscFunctionBegin; 1479 if (!*sp) PetscFunctionReturn(0); 1480 PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1); 1481 1482 if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 1483 ((PetscObject) (*sp))->refct = 0; 1484 /* if memory was published with AMS then destroy it */ 1485 ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr); 1486 1487 ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr); 1488 for (f = 0; f < dim; ++f) { 1489 /* ierr = PetscQuadratureDestroy((*sp)->functional[f]);CHKERRQ(ierr); */ 1490 } 1491 ierr = PetscFree((*sp)->functional);CHKERRQ(ierr); 1492 1493 ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr); 1494 ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 1495 PetscFunctionReturn(0); 1496 } 1497 1498 #undef __FUNCT__ 1499 #define __FUNCT__ "PetscDualSpaceCreate" 1500 /*@ 1501 PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType(). 1502 1503 Collective on MPI_Comm 1504 1505 Input Parameter: 1506 . comm - The communicator for the PetscDualSpace object 1507 1508 Output Parameter: 1509 . sp - The PetscDualSpace object 1510 1511 Level: beginner 1512 1513 .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE 1514 @*/ 1515 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp) 1516 { 1517 PetscDualSpace s; 1518 PetscErrorCode ierr; 1519 1520 PetscFunctionBegin; 1521 PetscValidPointer(sp, 2); 1522 *sp = NULL; 1523 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1524 ierr = PetscFEInitializePackage();CHKERRQ(ierr); 1525 #endif 1526 1527 ierr = PetscHeaderCreate(s, _p_PetscDualSpace, struct _PetscDualSpaceOps, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr); 1528 ierr = PetscMemzero(s->ops, sizeof(struct _PetscDualSpaceOps));CHKERRQ(ierr); 1529 1530 s->order = 0; 1531 1532 *sp = s; 1533 PetscFunctionReturn(0); 1534 } 1535 1536 #undef __FUNCT__ 1537 #define __FUNCT__ "PetscDualSpaceGetDM" 1538 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm) 1539 { 1540 PetscFunctionBegin; 1541 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1542 PetscValidPointer(dm, 2); 1543 *dm = sp->dm; 1544 PetscFunctionReturn(0); 1545 } 1546 1547 #undef __FUNCT__ 1548 #define __FUNCT__ "PetscDualSpaceSetDM" 1549 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm) 1550 { 1551 PetscFunctionBegin; 1552 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1553 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1554 sp->dm = dm; 1555 PetscFunctionReturn(0); 1556 } 1557 1558 #undef __FUNCT__ 1559 #define __FUNCT__ "PetscDualSpaceGetOrder" 1560 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order) 1561 { 1562 PetscFunctionBegin; 1563 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1564 PetscValidPointer(order, 2); 1565 *order = sp->order; 1566 PetscFunctionReturn(0); 1567 } 1568 1569 #undef __FUNCT__ 1570 #define __FUNCT__ "PetscDualSpaceSetOrder" 1571 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order) 1572 { 1573 PetscFunctionBegin; 1574 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1575 sp->order = order; 1576 PetscFunctionReturn(0); 1577 } 1578 1579 #undef __FUNCT__ 1580 #define __FUNCT__ "PetscDualSpaceGetFunctional" 1581 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional) 1582 { 1583 PetscInt dim; 1584 PetscErrorCode ierr; 1585 1586 PetscFunctionBegin; 1587 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1588 PetscValidPointer(functional, 3); 1589 ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 1590 if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim); 1591 *functional = sp->functional[i]; 1592 PetscFunctionReturn(0); 1593 } 1594 1595 #undef __FUNCT__ 1596 #define __FUNCT__ "PetscDualSpaceGetDimension" 1597 /* Dimension of the space, i.e. number of basis vectors */ 1598 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim) 1599 { 1600 PetscErrorCode ierr; 1601 1602 PetscFunctionBegin; 1603 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1604 PetscValidPointer(dim, 2); 1605 *dim = 0; 1606 if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 1607 PetscFunctionReturn(0); 1608 } 1609 1610 #undef __FUNCT__ 1611 #define __FUNCT__ "PetscDualSpaceSetUp_Lagrange" 1612 PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 1613 { 1614 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 1615 DM dm = sp->dm; 1616 PetscInt order = sp->order; 1617 PetscSection csection; 1618 Vec coordinates; 1619 PetscScalar *coords; 1620 PetscInt pdim, dim, vStart, vEnd, cStart, coneSize, f = 0; 1621 PetscErrorCode ierr; 1622 1623 PetscFunctionBegin; 1624 ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 1625 ierr = PetscMalloc(pdim * sizeof(PetscQuadrature), &sp->functional);CHKERRQ(ierr); 1626 /* Classify element type */ 1627 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1628 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 1629 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1630 ierr = DMPlexGetConeSize(dm, cStart, &coneSize);CHKERRQ(ierr); 1631 ierr = DMPlexGetCoordinateSection(dm, &csection);CHKERRQ(ierr); 1632 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1633 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1634 if (coneSize == dim+1) { 1635 PetscInt *closure = NULL, closureSize, c; 1636 1637 /* Simplex */ 1638 if (order > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle order %d basis", order); 1639 ierr = DMPlexGetTransitiveClosure(dm, cStart, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1640 for (c = 0; c < closureSize*2; c += 2) { 1641 const PetscInt p = closure[c]; 1642 /* Vertices */ 1643 if ((p >= vStart) && (p < vEnd)) { 1644 PetscReal *qpoints, *qweights; 1645 PetscInt dof, off, d; 1646 1647 sp->functional[f].numQuadPoints = 1; 1648 ierr = PetscMalloc(sp->functional[f].numQuadPoints*dim * sizeof(PetscReal), &qpoints);CHKERRQ(ierr); 1649 ierr = PetscMalloc(sp->functional[f].numQuadPoints * sizeof(PetscReal), &qweights);CHKERRQ(ierr); 1650 ierr = PetscSectionGetDof(csection, p, &dof);CHKERRQ(ierr); 1651 ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr); 1652 if (dof != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of coordinates %d does not match spatial dimension %d", dof, dim); 1653 for (d = 0; d < dof; ++d) {qpoints[d] = coords[off+d];} 1654 qweights[0] = 1.0; 1655 sp->functional[f].quadPoints = qpoints; 1656 sp->functional[f].quadWeights = qweights; 1657 ++f; 1658 } 1659 } 1660 ierr = DMPlexRestoreTransitiveClosure(dm, cStart, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1661 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle cells with cone size %d", coneSize); 1662 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1663 PetscFunctionReturn(0); 1664 } 1665 1666 #undef __FUNCT__ 1667 #define __FUNCT__ "PetscDualSpaceGetDimension_Lagrange" 1668 PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim) 1669 { 1670 PetscInt deg = sp->order; 1671 PetscReal D = 1.0; 1672 PetscInt n, i; 1673 PetscErrorCode ierr; 1674 1675 PetscFunctionBegin; 1676 /* TODO: Assumes simplices */ 1677 ierr = DMPlexGetDimension(sp->dm, &n);CHKERRQ(ierr); 1678 for (i = 1; i <= n; ++i) { 1679 D *= ((PetscReal) (deg+i))/i; 1680 } 1681 *dim = (PetscInt) (D + 0.5); 1682 PetscFunctionReturn(0); 1683 } 1684 1685 #undef __FUNCT__ 1686 #define __FUNCT__ "PetscDualSpaceInitialize_Lagrange" 1687 PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 1688 { 1689 PetscFunctionBegin; 1690 sp->ops->setfromoptions = NULL; 1691 sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 1692 sp->ops->view = NULL; 1693 sp->ops->destroy = NULL; 1694 sp->ops->getdimension = PetscDualSpaceGetDimension_Lagrange; 1695 PetscFunctionReturn(0); 1696 } 1697 1698 /*MC 1699 PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals 1700 1701 Level: intermediate 1702 1703 .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 1704 M*/ 1705 1706 #undef __FUNCT__ 1707 #define __FUNCT__ "PetscDualSpaceCreate_Lagrange" 1708 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 1709 { 1710 PetscDualSpace_Lag *lag; 1711 PetscErrorCode ierr; 1712 1713 PetscFunctionBegin; 1714 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1715 ierr = PetscNewLog(sp, PetscDualSpace_Lag, &lag);CHKERRQ(ierr); 1716 sp->data = lag; 1717 1718 /* lag->n = 0; */ 1719 1720 ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr); 1721 PetscFunctionReturn(0); 1722 } 1723 1724 1725 PetscInt PETSCFE_CLASSID = 0; 1726 1727 #undef __FUNCT__ 1728 #define __FUNCT__ "PetscFEView" 1729 /*@C 1730 PetscFEView - Views a PetscFE 1731 1732 Collective on PetscFE 1733 1734 Input Parameter: 1735 + sp - the PetscFE object to view 1736 - v - the viewer 1737 1738 Level: developer 1739 1740 .seealso PetscFEDestroy() 1741 @*/ 1742 PetscErrorCode PetscFEView(PetscFE fem, PetscViewer v) 1743 { 1744 PetscErrorCode ierr; 1745 1746 PetscFunctionBegin; 1747 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1748 if (!v) { 1749 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &v);CHKERRQ(ierr); 1750 } 1751 if (fem->ops->view) { 1752 ierr = (*fem->ops->view)(fem, v);CHKERRQ(ierr); 1753 } 1754 PetscFunctionReturn(0); 1755 } 1756 1757 #undef __FUNCT__ 1758 #define __FUNCT__ "PetscFEDestroy" 1759 /*@ 1760 PetscFEDestroy - Destroys a PetscFE object 1761 1762 Collective on PetscFE 1763 1764 Input Parameter: 1765 . fem - the PetscFE object to destroy 1766 1767 Level: developer 1768 1769 .seealso PetscFEView() 1770 @*/ 1771 PetscErrorCode PetscFEDestroy(PetscFE *fem) 1772 { 1773 PetscErrorCode ierr; 1774 1775 PetscFunctionBegin; 1776 if (!*fem) PetscFunctionReturn(0); 1777 PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1); 1778 1779 if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);} 1780 ((PetscObject) (*fem))->refct = 0; 1781 /* if memory was published with AMS then destroy it */ 1782 ierr = PetscObjectAMSViewOff((PetscObject) *fem);CHKERRQ(ierr); 1783 1784 ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr); 1785 ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr); 1786 PetscFunctionReturn(0); 1787 } 1788 1789 #undef __FUNCT__ 1790 #define __FUNCT__ "PetscFECreate" 1791 /*@ 1792 PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType(). 1793 1794 Collective on MPI_Comm 1795 1796 Input Parameter: 1797 . comm - The communicator for the PetscFE object 1798 1799 Output Parameter: 1800 . fem - The PetscFE object 1801 1802 Level: beginner 1803 1804 .seealso: PetscFESetType(), PETSCFEGALERKIN 1805 @*/ 1806 PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem) 1807 { 1808 PetscFE f; 1809 PetscErrorCode ierr; 1810 1811 PetscFunctionBegin; 1812 PetscValidPointer(fem, 2); 1813 *fem = NULL; 1814 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1815 ierr = PetscFEInitializePackage();CHKERRQ(ierr); 1816 #endif 1817 1818 ierr = PetscHeaderCreate(f, _p_PetscFE, struct _PetscFEOps, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr); 1819 ierr = PetscMemzero(f->ops, sizeof(struct _PetscFEOps));CHKERRQ(ierr); 1820 1821 f->basisSpace = NULL; 1822 f->dualSpace = NULL; 1823 1824 *fem = f; 1825 PetscFunctionReturn(0); 1826 } 1827 1828 #undef __FUNCT__ 1829 #define __FUNCT__ "PetscFEGetDimension" 1830 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 1831 { 1832 PetscErrorCode ierr; 1833 1834 PetscFunctionBegin; 1835 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1836 PetscValidPointer(dim, 2); 1837 ierr = PetscSpaceGetDimension(fem->basisSpace, dim);CHKERRQ(ierr); 1838 PetscFunctionReturn(0); 1839 } 1840 1841 #undef __FUNCT__ 1842 #define __FUNCT__ "PetscFEGetNumComponents" 1843 PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp) 1844 { 1845 PetscFunctionBegin; 1846 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1847 PetscValidPointer(comp, 2); 1848 *comp = 1; 1849 PetscFunctionReturn(0); 1850 } 1851 1852 #undef __FUNCT__ 1853 #define __FUNCT__ "PetscFEGetBasisSpace" 1854 PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp) 1855 { 1856 PetscFunctionBegin; 1857 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1858 PetscValidPointer(sp, 2); 1859 *sp = fem->basisSpace; 1860 PetscFunctionReturn(0); 1861 } 1862 1863 #undef __FUNCT__ 1864 #define __FUNCT__ "PetscFESetBasisSpace" 1865 PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp) 1866 { 1867 PetscFunctionBegin; 1868 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1869 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2); 1870 fem->basisSpace = sp; 1871 PetscFunctionReturn(0); 1872 } 1873 1874 #undef __FUNCT__ 1875 #define __FUNCT__ "PetscFEGetDualSpace" 1876 PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp) 1877 { 1878 PetscFunctionBegin; 1879 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1880 PetscValidPointer(sp, 2); 1881 *sp = fem->dualSpace; 1882 PetscFunctionReturn(0); 1883 } 1884 1885 #undef __FUNCT__ 1886 #define __FUNCT__ "PetscFESetDualSpace" 1887 PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp) 1888 { 1889 PetscFunctionBegin; 1890 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1891 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 1892 fem->dualSpace = sp; 1893 PetscFunctionReturn(0); 1894 } 1895 1896 #undef __FUNCT__ 1897 #define __FUNCT__ "PetscFEGetTabulation" 1898 PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 1899 { 1900 DM dm; 1901 PetscReal *tmpB, *invV; 1902 PetscInt pdim; /* Dimension of FE space P */ 1903 PetscInt dim; /* Spatial dimension */ 1904 PetscInt p, j, k; 1905 PetscErrorCode ierr; 1906 1907 PetscFunctionBegin; 1908 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1909 PetscValidPointer(points, 3); 1910 if (B) PetscValidPointer(B, 4); 1911 if (D) PetscValidPointer(D, 5); 1912 if (H) PetscValidPointer(H, 6); 1913 ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 1914 1915 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1916 ierr = PetscSpaceGetDimension(fem->basisSpace, &pdim);CHKERRQ(ierr); 1917 /* if (nvalues%dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate values %d must be divisible by the spatial dimension %d", nvalues, dim); */ 1918 1919 if (pdim != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Code only works for linear elements right now"); 1920 1921 if (B) { 1922 ierr = DMGetWorkArray(dm, npoints*pdim, PETSC_REAL, B);CHKERRQ(ierr); 1923 ierr = DMGetWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr); 1924 } 1925 if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*dim, PETSC_REAL, D);CHKERRQ(ierr);} 1926 if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*dim*dim, PETSC_REAL, H);CHKERRQ(ierr);} 1927 ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr); 1928 1929 ierr = DMGetWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr); 1930 for (j = 0; j < pdim; ++j) { 1931 PetscReal *Bf; 1932 PetscQuadrature f; 1933 PetscInt q; 1934 1935 ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f); 1936 ierr = DMGetWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr); 1937 ierr = PetscSpaceEvaluate(fem->basisSpace, f.numQuadPoints, f.quadPoints, Bf, NULL, NULL);CHKERRQ(ierr); 1938 for (k = 0; k < pdim; ++k) { 1939 /* n_j \cdot \phi_k */ 1940 invV[j*pdim+k] = 0.0; 1941 for (q = 0; q < f.numQuadPoints; ++q) { 1942 invV[j*pdim+k] += Bf[q*pdim+k]*f.quadWeights[q]; 1943 } 1944 } 1945 ierr = DMRestoreWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr); 1946 } 1947 { 1948 PetscReal *work; 1949 PetscBLASInt *pivots; 1950 PetscBLASInt n = pdim, info; 1951 1952 ierr = DMGetWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr); 1953 ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr); 1954 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invV, &n, pivots, &info)); 1955 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invV, &n, pivots, work, &n, &info)); 1956 ierr = DMRestoreWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr); 1957 ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr); 1958 } 1959 for (p = 0; p < npoints; ++p) { 1960 if (B) { 1961 /* Multiply by V^{-1} (pdim x pdim) */ 1962 for (j = 0; j < pdim; ++j) { 1963 const PetscInt i = p*pdim + j; 1964 1965 (*B)[i] = 0.0; 1966 for (k = 0; k < pdim; ++k) { 1967 (*B)[i] += invV[k*pdim+j] * tmpB[p*pdim + k]; 1968 } 1969 } 1970 } 1971 } 1972 ierr = DMRestoreWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr); 1973 if (B) { 1974 ierr = DMRestoreWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr); 1975 } 1976 PetscFunctionReturn(0); 1977 } 1978 1979 #undef __FUNCT__ 1980 #define __FUNCT__ "PetscFERestoreTabulation" 1981 PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **D2) 1982 { 1983 DM dm; 1984 PetscErrorCode ierr; 1985 1986 PetscFunctionBegin; 1987 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1988 ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 1989 if (B) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, B);CHKERRQ(ierr);} 1990 if (D) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D);CHKERRQ(ierr);} 1991 if (D2) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D2);CHKERRQ(ierr);} 1992 PetscFunctionReturn(0); 1993 } 1994