137045ce4SJed Brown /* Discretization tools */ 237045ce4SJed Brown 3a6fc04d9SSatish Balay #include <petscconf.h> 4a6fc04d9SSatish Balay #if defined(PETSC_HAVE_MATHIMF_H) 5a6fc04d9SSatish Balay #include <mathimf.h> /* this needs to be included before math.h */ 6a6fc04d9SSatish Balay #endif 7a6fc04d9SSatish Balay 837045ce4SJed Brown #include <petscdt.h> /*I "petscdt.h" I*/ 937045ce4SJed Brown #include <petscblaslapack.h> 10194825f6SJed Brown #include <petsc-private/petscimpl.h> 11665c2dedSJed Brown #include <petscviewer.h> 1237045ce4SJed Brown 1337045ce4SJed Brown #undef __FUNCT__ 1437045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval" 1537045ce4SJed Brown /*@ 1637045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 1737045ce4SJed Brown 1837045ce4SJed Brown Not Collective 1937045ce4SJed Brown 2037045ce4SJed Brown Input Arguments: 2137045ce4SJed Brown + npoints - number of spatial points to evaluate at 2237045ce4SJed Brown . points - array of locations to evaluate at 2337045ce4SJed Brown . ndegree - number of basis degrees to evaluate 2437045ce4SJed Brown - degrees - sorted array of degrees to evaluate 2537045ce4SJed Brown 2637045ce4SJed Brown Output Arguments: 270298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 280298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 290298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 3037045ce4SJed Brown 3137045ce4SJed Brown Level: intermediate 3237045ce4SJed Brown 3337045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 3437045ce4SJed Brown @*/ 3537045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 3637045ce4SJed Brown { 3737045ce4SJed Brown PetscInt i,maxdegree; 3837045ce4SJed Brown 3937045ce4SJed Brown PetscFunctionBegin; 4037045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 4137045ce4SJed Brown maxdegree = degrees[ndegree-1]; 4237045ce4SJed Brown for (i=0; i<npoints; i++) { 4337045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 4437045ce4SJed Brown PetscInt j,k; 4537045ce4SJed Brown x = points[i]; 4637045ce4SJed Brown pm2 = 0; 4737045ce4SJed Brown pm1 = 1; 4837045ce4SJed Brown pd2 = 0; 4937045ce4SJed Brown pd1 = 0; 5037045ce4SJed Brown pdd2 = 0; 5137045ce4SJed Brown pdd1 = 0; 5237045ce4SJed Brown k = 0; 5337045ce4SJed Brown if (degrees[k] == 0) { 5437045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 5537045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 5637045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 5737045ce4SJed Brown k++; 5837045ce4SJed Brown } 5937045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 6037045ce4SJed Brown PetscReal p,d,dd; 6137045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 6237045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 6337045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 6437045ce4SJed Brown pm2 = pm1; 6537045ce4SJed Brown pm1 = p; 6637045ce4SJed Brown pd2 = pd1; 6737045ce4SJed Brown pd1 = d; 6837045ce4SJed Brown pdd2 = pdd1; 6937045ce4SJed Brown pdd1 = dd; 7037045ce4SJed Brown if (degrees[k] == j) { 7137045ce4SJed Brown if (B) B[i*ndegree+k] = p; 7237045ce4SJed Brown if (D) D[i*ndegree+k] = d; 7337045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 7437045ce4SJed Brown } 7537045ce4SJed Brown } 7637045ce4SJed Brown } 7737045ce4SJed Brown PetscFunctionReturn(0); 7837045ce4SJed Brown } 7937045ce4SJed Brown 8037045ce4SJed Brown #undef __FUNCT__ 8137045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature" 8237045ce4SJed Brown /*@ 8337045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 8437045ce4SJed Brown 8537045ce4SJed Brown Not Collective 8637045ce4SJed Brown 8737045ce4SJed Brown Input Arguments: 8837045ce4SJed Brown + npoints - number of points 8937045ce4SJed Brown . a - left end of interval (often-1) 9037045ce4SJed Brown - b - right end of interval (often +1) 9137045ce4SJed Brown 9237045ce4SJed Brown Output Arguments: 9337045ce4SJed Brown + x - quadrature points 9437045ce4SJed Brown - w - quadrature weights 9537045ce4SJed Brown 9637045ce4SJed Brown Level: intermediate 9737045ce4SJed Brown 9837045ce4SJed Brown References: 9937045ce4SJed Brown Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 10037045ce4SJed Brown 10137045ce4SJed Brown .seealso: PetscDTLegendreEval() 10237045ce4SJed Brown @*/ 10337045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 10437045ce4SJed Brown { 10537045ce4SJed Brown PetscErrorCode ierr; 10637045ce4SJed Brown PetscInt i; 10737045ce4SJed Brown PetscReal *work; 10837045ce4SJed Brown PetscScalar *Z; 10937045ce4SJed Brown PetscBLASInt N,LDZ,info; 11037045ce4SJed Brown 11137045ce4SJed Brown PetscFunctionBegin; 11237045ce4SJed Brown /* Set up the Golub-Welsch system */ 11337045ce4SJed Brown for (i=0; i<npoints; i++) { 11437045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 11537045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 11637045ce4SJed Brown } 11737045ce4SJed Brown ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 11837045ce4SJed Brown ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr); 119c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 12037045ce4SJed Brown LDZ = N; 12137045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1228b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 12337045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 1241c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 12537045ce4SJed Brown 12637045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 12737045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 12837045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 12937045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 1300d644c17SKarl Rupp 13137045ce4SJed Brown w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints])); 13237045ce4SJed Brown } 13337045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 13437045ce4SJed Brown PetscFunctionReturn(0); 13537045ce4SJed Brown } 136194825f6SJed Brown 137194825f6SJed Brown #undef __FUNCT__ 138494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal" 139494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 140494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 141494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 142494e7359SMatthew G. Knepley { 143494e7359SMatthew G. Knepley PetscReal f = 1.0; 144494e7359SMatthew G. Knepley PetscInt i; 145494e7359SMatthew G. Knepley 146494e7359SMatthew G. Knepley PetscFunctionBegin; 147494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 148494e7359SMatthew G. Knepley *factorial = f; 149494e7359SMatthew G. Knepley PetscFunctionReturn(0); 150494e7359SMatthew G. Knepley } 151494e7359SMatthew G. Knepley 152494e7359SMatthew G. Knepley #undef __FUNCT__ 153494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi" 154494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 155494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 156494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 157494e7359SMatthew G. Knepley { 158494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 159494e7359SMatthew G. Knepley PetscInt k; 160494e7359SMatthew G. Knepley 161494e7359SMatthew G. Knepley PetscFunctionBegin; 162494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 163494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 164494e7359SMatthew G. Knepley apb = a + b; 165494e7359SMatthew G. Knepley pn2 = 1.0; 166494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 167494e7359SMatthew G. Knepley *P = 0.0; 168494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 169494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 170494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 171494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 172494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 173494e7359SMatthew G. Knepley 174494e7359SMatthew G. Knepley a2 = a2 / a1; 175494e7359SMatthew G. Knepley a3 = a3 / a1; 176494e7359SMatthew G. Knepley a4 = a4 / a1; 177494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 178494e7359SMatthew G. Knepley pn2 = pn1; 179494e7359SMatthew G. Knepley pn1 = *P; 180494e7359SMatthew G. Knepley } 181494e7359SMatthew G. Knepley PetscFunctionReturn(0); 182494e7359SMatthew G. Knepley } 183494e7359SMatthew G. Knepley 184494e7359SMatthew G. Knepley #undef __FUNCT__ 185494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative" 186494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 187494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 188494e7359SMatthew G. Knepley { 189494e7359SMatthew G. Knepley PetscReal nP; 190494e7359SMatthew G. Knepley PetscErrorCode ierr; 191494e7359SMatthew G. Knepley 192494e7359SMatthew G. Knepley PetscFunctionBegin; 193494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 194494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 195494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 196494e7359SMatthew G. Knepley PetscFunctionReturn(0); 197494e7359SMatthew G. Knepley } 198494e7359SMatthew G. Knepley 199494e7359SMatthew G. Knepley #undef __FUNCT__ 200494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 201494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 202494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 203494e7359SMatthew G. Knepley { 204494e7359SMatthew G. Knepley PetscFunctionBegin; 205494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 206494e7359SMatthew G. Knepley *eta = y; 207494e7359SMatthew G. Knepley PetscFunctionReturn(0); 208494e7359SMatthew G. Knepley } 209494e7359SMatthew G. Knepley 210494e7359SMatthew G. Knepley #undef __FUNCT__ 211494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 212494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 213494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 214494e7359SMatthew G. Knepley { 215494e7359SMatthew G. Knepley PetscFunctionBegin; 216494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 217494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 218494e7359SMatthew G. Knepley *zeta = z; 219494e7359SMatthew G. Knepley PetscFunctionReturn(0); 220494e7359SMatthew G. Knepley } 221494e7359SMatthew G. Knepley 222494e7359SMatthew G. Knepley #undef __FUNCT__ 223494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 224494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 225494e7359SMatthew G. Knepley { 226494e7359SMatthew G. Knepley PetscInt maxIter = 100; 227494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 228a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 229494e7359SMatthew G. Knepley PetscInt k; 230494e7359SMatthew G. Knepley PetscErrorCode ierr; 231494e7359SMatthew G. Knepley 232494e7359SMatthew G. Knepley PetscFunctionBegin; 233a8291ba1SSatish Balay 234a8291ba1SSatish Balay a1 = pow(2, a+b+1); 235a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 236a8291ba1SSatish Balay a2 = tgamma(a + npoints + 1); 237a8291ba1SSatish Balay a3 = tgamma(b + npoints + 1); 238a8291ba1SSatish Balay a4 = tgamma(a + b + npoints + 1); 239a8291ba1SSatish Balay #else 240a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 241a8291ba1SSatish Balay #endif 242a8291ba1SSatish Balay 243494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 244494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 245494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 246494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 247494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 2487f1c68b3SMatthew G. Knepley PetscReal r = -cos((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 249494e7359SMatthew G. Knepley PetscInt j; 250494e7359SMatthew G. Knepley 251494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 252494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 253494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 254494e7359SMatthew G. Knepley PetscInt i; 255494e7359SMatthew G. Knepley 256494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 257494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 258494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 259494e7359SMatthew G. Knepley delta = f / (fp - f * s); 260494e7359SMatthew G. Knepley r = r - delta; 261494e7359SMatthew G. Knepley if (fabs(delta) < eps) break; 262494e7359SMatthew G. Knepley } 263494e7359SMatthew G. Knepley x[k] = r; 264494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 265494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 266494e7359SMatthew G. Knepley } 267494e7359SMatthew G. Knepley PetscFunctionReturn(0); 268494e7359SMatthew G. Knepley } 269494e7359SMatthew G. Knepley 270494e7359SMatthew G. Knepley #undef __FUNCT__ 271494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 272fd9d31fbSMatthew G. Knepley /*@C 273494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 274494e7359SMatthew G. Knepley 275494e7359SMatthew G. Knepley Not Collective 276494e7359SMatthew G. Knepley 277494e7359SMatthew G. Knepley Input Arguments: 278494e7359SMatthew G. Knepley + dim - The simplex dimension 279494e7359SMatthew G. Knepley . npoints - number of points 280494e7359SMatthew G. Knepley . a - left end of interval (often-1) 281494e7359SMatthew G. Knepley - b - right end of interval (often +1) 282494e7359SMatthew G. Knepley 283494e7359SMatthew G. Knepley Output Arguments: 284494e7359SMatthew G. Knepley + points - quadrature points 285494e7359SMatthew G. Knepley - weights - quadrature weights 286494e7359SMatthew G. Knepley 287494e7359SMatthew G. Knepley Level: intermediate 288494e7359SMatthew G. Knepley 289494e7359SMatthew G. Knepley References: 290494e7359SMatthew G. Knepley Karniadakis and Sherwin. 291494e7359SMatthew G. Knepley FIAT 292494e7359SMatthew G. Knepley 293494e7359SMatthew G. Knepley .seealso: PetscDTGaussQuadrature() 294494e7359SMatthew G. Knepley @*/ 295494e7359SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscReal *points[], PetscReal *weights[]) 296494e7359SMatthew G. Knepley { 297494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 298494e7359SMatthew G. Knepley PetscInt i, j, k; 299494e7359SMatthew G. Knepley PetscErrorCode ierr; 300494e7359SMatthew G. Knepley 301494e7359SMatthew G. Knepley PetscFunctionBegin; 302494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 303494e7359SMatthew G. Knepley switch (dim) { 304494e7359SMatthew G. Knepley case 1: 305494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints * sizeof(PetscReal), &x);CHKERRQ(ierr); 306494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 307494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, w);CHKERRQ(ierr); 308494e7359SMatthew G. Knepley break; 309494e7359SMatthew G. Knepley case 2: 310494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints*npoints*2 * sizeof(PetscReal), &x);CHKERRQ(ierr); 311494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints*npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 312494e7359SMatthew G. Knepley ierr = PetscMalloc4(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy);CHKERRQ(ierr); 313494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 314494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 315494e7359SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 316494e7359SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 317494e7359SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 318494e7359SMatthew G. Knepley w[i*npoints+j] = 0.5 * wx[i] * wy[j]; 319494e7359SMatthew G. Knepley } 320494e7359SMatthew G. Knepley } 321494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 322494e7359SMatthew G. Knepley break; 323494e7359SMatthew G. Knepley case 3: 324494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints*npoints*3 * sizeof(PetscReal), &x);CHKERRQ(ierr); 325494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints*npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 326494e7359SMatthew G. Knepley ierr = PetscMalloc6(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy,npoints,PetscReal,&pz,npoints,PetscReal,&wz);CHKERRQ(ierr); 327494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 328494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 329494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 330494e7359SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 331494e7359SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 332494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 333494e7359SMatthew G. Knepley 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); 334494e7359SMatthew G. Knepley w[(i*npoints+j)*npoints+k] = 0.125 * wx[i] * wy[j] * wz[k]; 335494e7359SMatthew G. Knepley } 336494e7359SMatthew G. Knepley } 337494e7359SMatthew G. Knepley } 338494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 339494e7359SMatthew G. Knepley break; 340494e7359SMatthew G. Knepley default: 341494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 342494e7359SMatthew G. Knepley } 343494e7359SMatthew G. Knepley if (points) *points = x; 344494e7359SMatthew G. Knepley if (weights) *weights = w; 345494e7359SMatthew G. Knepley PetscFunctionReturn(0); 346494e7359SMatthew G. Knepley } 347494e7359SMatthew G. Knepley 348494e7359SMatthew G. Knepley #undef __FUNCT__ 349194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR" 350194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 351194825f6SJed Brown * A in column-major format 352194825f6SJed Brown * Ainv in row-major format 353194825f6SJed Brown * tau has length m 354194825f6SJed Brown * worksize must be >= max(1,n) 355194825f6SJed Brown */ 356194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 357194825f6SJed Brown { 358194825f6SJed Brown PetscErrorCode ierr; 359194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 360194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 361194825f6SJed Brown 362194825f6SJed Brown PetscFunctionBegin; 363194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 364194825f6SJed Brown { 365194825f6SJed Brown PetscInt i,j; 366194825f6SJed Brown ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr); 367194825f6SJed Brown for (j=0; j<n; j++) { 368194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 369194825f6SJed Brown } 370194825f6SJed Brown mstride = m; 371194825f6SJed Brown } 372194825f6SJed Brown #else 373194825f6SJed Brown A = A_in; 374194825f6SJed Brown Ainv = Ainv_out; 375194825f6SJed Brown #endif 376194825f6SJed Brown 377194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 378194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 379194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 380194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 381194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 382194825f6SJed Brown LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info); 383194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 384194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 385194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 386194825f6SJed Brown 387194825f6SJed Brown /* Extract an explicit representation of Q */ 388194825f6SJed Brown Q = Ainv; 389194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 390194825f6SJed Brown K = N; /* full rank */ 391194825f6SJed Brown LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info); 392194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 393194825f6SJed Brown 394194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 395194825f6SJed Brown Alpha = 1.0; 396194825f6SJed Brown ldb = lda; 397194825f6SJed Brown BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb); 398194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 399194825f6SJed Brown 400194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 401194825f6SJed Brown { 402194825f6SJed Brown PetscInt i; 403194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 404194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 405194825f6SJed Brown } 406194825f6SJed Brown #endif 407194825f6SJed Brown PetscFunctionReturn(0); 408194825f6SJed Brown } 409194825f6SJed Brown 410194825f6SJed Brown #undef __FUNCT__ 411194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate" 412194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 413194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 414194825f6SJed Brown { 415194825f6SJed Brown PetscErrorCode ierr; 416194825f6SJed Brown PetscReal *Bv; 417194825f6SJed Brown PetscInt i,j; 418194825f6SJed Brown 419194825f6SJed Brown PetscFunctionBegin; 420194825f6SJed Brown ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr); 421194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 422194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 423194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 424194825f6SJed Brown for (i=0; i<ninterval; i++) { 425194825f6SJed Brown for (j=0; j<ndegree; j++) { 426194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 427194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 428194825f6SJed Brown } 429194825f6SJed Brown } 430194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 431194825f6SJed Brown PetscFunctionReturn(0); 432194825f6SJed Brown } 433194825f6SJed Brown 434194825f6SJed Brown #undef __FUNCT__ 435194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly" 436194825f6SJed Brown /*@ 437194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 438194825f6SJed Brown 439194825f6SJed Brown Not Collective 440194825f6SJed Brown 441194825f6SJed Brown Input Arguments: 442194825f6SJed Brown + degree - degree of reconstruction polynomial 443194825f6SJed Brown . nsource - number of source intervals 444194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 445194825f6SJed Brown . ntarget - number of target intervals 446194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 447194825f6SJed Brown 448194825f6SJed Brown Output Arguments: 449194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 450194825f6SJed Brown 451194825f6SJed Brown Level: advanced 452194825f6SJed Brown 453194825f6SJed Brown .seealso: PetscDTLegendreEval() 454194825f6SJed Brown @*/ 455194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 456194825f6SJed Brown { 457194825f6SJed Brown PetscErrorCode ierr; 458194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 459194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 460194825f6SJed Brown PetscScalar *tau,*work; 461194825f6SJed Brown 462194825f6SJed Brown PetscFunctionBegin; 463194825f6SJed Brown PetscValidRealPointer(sourcex,3); 464194825f6SJed Brown PetscValidRealPointer(targetx,5); 465194825f6SJed Brown PetscValidRealPointer(R,6); 466194825f6SJed Brown 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); 467194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 468194825f6SJed Brown for (i=0; i<nsource; i++) { 469194825f6SJed Brown 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]); 470194825f6SJed Brown } 471194825f6SJed Brown for (i=0; i<ntarget; i++) { 472194825f6SJed Brown 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]); 473194825f6SJed Brown } 474194825f6SJed Brown #endif 475194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 476194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 477194825f6SJed Brown center = (xmin + xmax)/2; 478194825f6SJed Brown hscale = (xmax - xmin)/2; 479194825f6SJed Brown worksize = nsource; 480194825f6SJed Brown ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr); 48182772646SJed Brown ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr); 482194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 483194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 484194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 485194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 486194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 487194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 488194825f6SJed Brown for (i=0; i<ntarget; i++) { 489194825f6SJed Brown PetscReal rowsum = 0; 490194825f6SJed Brown for (j=0; j<nsource; j++) { 491194825f6SJed Brown PetscReal sum = 0; 492194825f6SJed Brown for (k=0; k<degree+1; k++) { 493194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 494194825f6SJed Brown } 495194825f6SJed Brown R[i*nsource+j] = sum; 496194825f6SJed Brown rowsum += sum; 497194825f6SJed Brown } 498194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 499194825f6SJed Brown } 500194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 501194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 502194825f6SJed Brown PetscFunctionReturn(0); 503194825f6SJed Brown } 504*29d8e7bcSMatthew G. Knepley 505*29d8e7bcSMatthew G. Knepley /* Basis Jet Tabulation 506*29d8e7bcSMatthew G. Knepley 507*29d8e7bcSMatthew G. Knepley We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We 508*29d8e7bcSMatthew G. Knepley follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can 509*29d8e7bcSMatthew G. Knepley be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis 510*29d8e7bcSMatthew G. Knepley as a prime basis. 511*29d8e7bcSMatthew G. Knepley 512*29d8e7bcSMatthew G. Knepley \psi_i = \sum_k \alpha_{ki} \phi_k 513*29d8e7bcSMatthew G. Knepley 514*29d8e7bcSMatthew G. Knepley Our nodal basis is defined in terms of the dual basis $n_j$ 515*29d8e7bcSMatthew G. Knepley 516*29d8e7bcSMatthew G. Knepley n_j \cdot \psi_i = \delta_{ji} 517*29d8e7bcSMatthew G. Knepley 518*29d8e7bcSMatthew G. Knepley and we may act on the first equation to obtain 519*29d8e7bcSMatthew G. Knepley 520*29d8e7bcSMatthew G. Knepley n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k 521*29d8e7bcSMatthew G. Knepley \delta_{ji} = \sum_k \alpha_{ki} V_{jk} 522*29d8e7bcSMatthew G. Knepley I = V \alpha 523*29d8e7bcSMatthew G. Knepley 524*29d8e7bcSMatthew G. Knepley so the coefficients of the nodal basis in the prime basis are 525*29d8e7bcSMatthew G. Knepley 526*29d8e7bcSMatthew G. Knepley \alpha = V^{-1} 527*29d8e7bcSMatthew G. Knepley 528*29d8e7bcSMatthew G. Knepley We will define the dual basis vectors $n_j$ using a quadrature rule. 529*29d8e7bcSMatthew G. Knepley 530*29d8e7bcSMatthew G. Knepley Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials 531*29d8e7bcSMatthew G. Knepley (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can 532*29d8e7bcSMatthew G. Knepley be implemented exactly as in FIAT using functionals $L_j$. 533*29d8e7bcSMatthew G. Knepley 534*29d8e7bcSMatthew G. Knepley I will have to count the degrees correctly for the Legendre product when we are on simplices. 535*29d8e7bcSMatthew G. Knepley 536*29d8e7bcSMatthew G. Knepley We will have three objects: 537*29d8e7bcSMatthew G. Knepley - Space, P: this just need point evaluation I think 538*29d8e7bcSMatthew G. Knepley - 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*29d8e7bcSMatthew G. Knepley - FEM: This keeps {P, P', Q} 540*29d8e7bcSMatthew G. Knepley */ 541*29d8e7bcSMatthew G. Knepley #include <petsc-private/petscfeimpl.h> 542*29d8e7bcSMatthew G. Knepley 543*29d8e7bcSMatthew G. Knepley PetscInt PETSCSPACE_CLASSID = 0; 544*29d8e7bcSMatthew G. Knepley 545*29d8e7bcSMatthew G. Knepley PetscFunctionList PetscSpaceList = NULL; 546*29d8e7bcSMatthew G. Knepley PetscBool PetscSpaceRegisterAllCalled = PETSC_FALSE; 547*29d8e7bcSMatthew G. Knepley 548*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 549*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceRegister" 550*29d8e7bcSMatthew G. Knepley /*@C 551*29d8e7bcSMatthew G. Knepley PetscSpaceRegister - Adds a new PetscSpace implementation 552*29d8e7bcSMatthew G. Knepley 553*29d8e7bcSMatthew G. Knepley Not Collective 554*29d8e7bcSMatthew G. Knepley 555*29d8e7bcSMatthew G. Knepley Input Parameters: 556*29d8e7bcSMatthew G. Knepley + name - The name of a new user-defined creation routine 557*29d8e7bcSMatthew G. Knepley - create_func - The creation routine itself 558*29d8e7bcSMatthew G. Knepley 559*29d8e7bcSMatthew G. Knepley Notes: 560*29d8e7bcSMatthew G. Knepley PetscSpaceRegister() may be called multiple times to add several user-defined PetscSpaces 561*29d8e7bcSMatthew G. Knepley 562*29d8e7bcSMatthew G. Knepley Sample usage: 563*29d8e7bcSMatthew G. Knepley .vb 564*29d8e7bcSMatthew G. Knepley PetscSpaceRegister("my_space", MyPetscSpaceCreate); 565*29d8e7bcSMatthew G. Knepley .ve 566*29d8e7bcSMatthew G. Knepley 567*29d8e7bcSMatthew G. Knepley Then, your PetscSpace type can be chosen with the procedural interface via 568*29d8e7bcSMatthew G. Knepley .vb 569*29d8e7bcSMatthew G. Knepley PetscSpaceCreate(MPI_Comm, PetscSpace *); 570*29d8e7bcSMatthew G. Knepley PetscSpaceSetType(PetscSpace, "my_space"); 571*29d8e7bcSMatthew G. Knepley .ve 572*29d8e7bcSMatthew G. Knepley or at runtime via the option 573*29d8e7bcSMatthew G. Knepley .vb 574*29d8e7bcSMatthew G. Knepley -petscspace_type my_space 575*29d8e7bcSMatthew G. Knepley .ve 576*29d8e7bcSMatthew G. Knepley 577*29d8e7bcSMatthew G. Knepley Level: advanced 578*29d8e7bcSMatthew G. Knepley 579*29d8e7bcSMatthew G. Knepley .keywords: PetscSpace, register 580*29d8e7bcSMatthew G. Knepley .seealso: PetscSpaceRegisterAll(), PetscSpaceRegisterDestroy() 581*29d8e7bcSMatthew G. Knepley 582*29d8e7bcSMatthew G. Knepley @*/ 583*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscSpace)) 584*29d8e7bcSMatthew G. Knepley { 585*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 586*29d8e7bcSMatthew G. Knepley 587*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 588*29d8e7bcSMatthew G. Knepley ierr = PetscFunctionListAdd(&PetscSpaceList, sname, function);CHKERRQ(ierr); 589*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 590*29d8e7bcSMatthew G. Knepley } 591*29d8e7bcSMatthew G. Knepley 592*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 593*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetType" 594*29d8e7bcSMatthew G. Knepley /*@C 595*29d8e7bcSMatthew G. Knepley PetscSpaceSetType - Builds a particular PetscSpace 596*29d8e7bcSMatthew G. Knepley 597*29d8e7bcSMatthew G. Knepley Collective on PetscSpace 598*29d8e7bcSMatthew G. Knepley 599*29d8e7bcSMatthew G. Knepley Input Parameters: 600*29d8e7bcSMatthew G. Knepley + sp - The PetscSpace object 601*29d8e7bcSMatthew G. Knepley - name - The kind of space 602*29d8e7bcSMatthew G. Knepley 603*29d8e7bcSMatthew G. Knepley Options Database Key: 604*29d8e7bcSMatthew G. Knepley . -petscspace_type <type> - Sets the PetscSpace type; use -help for a list of available types 605*29d8e7bcSMatthew G. Knepley 606*29d8e7bcSMatthew G. Knepley Level: intermediate 607*29d8e7bcSMatthew G. Knepley 608*29d8e7bcSMatthew G. Knepley .keywords: PetscSpace, set, type 609*29d8e7bcSMatthew G. Knepley .seealso: PetscSpaceGetType(), PetscSpaceCreate() 610*29d8e7bcSMatthew G. Knepley @*/ 611*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetType(PetscSpace sp, PetscSpaceType name) 612*29d8e7bcSMatthew G. Knepley { 613*29d8e7bcSMatthew G. Knepley PetscErrorCode (*r)(PetscSpace); 614*29d8e7bcSMatthew G. Knepley PetscBool match; 615*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 616*29d8e7bcSMatthew G. Knepley 617*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 618*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 619*29d8e7bcSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 620*29d8e7bcSMatthew G. Knepley if (match) PetscFunctionReturn(0); 621*29d8e7bcSMatthew G. Knepley 622*29d8e7bcSMatthew G. Knepley if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 623*29d8e7bcSMatthew G. Knepley ierr = PetscFunctionListFind(PetscSpaceList, name, &r);CHKERRQ(ierr); 624*29d8e7bcSMatthew G. Knepley if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSpace type: %s", name); 625*29d8e7bcSMatthew G. Knepley 626*29d8e7bcSMatthew G. Knepley if (sp->ops->destroy) { 627*29d8e7bcSMatthew G. Knepley ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 628*29d8e7bcSMatthew G. Knepley sp->ops->destroy = NULL; 629*29d8e7bcSMatthew G. Knepley } 630*29d8e7bcSMatthew G. Knepley ierr = (*r)(sp);CHKERRQ(ierr); 631*29d8e7bcSMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 632*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 633*29d8e7bcSMatthew G. Knepley } 634*29d8e7bcSMatthew G. Knepley 635*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 636*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceGetType" 637*29d8e7bcSMatthew G. Knepley /*@C 638*29d8e7bcSMatthew G. Knepley PetscSpaceGetType - Gets the PetscSpace type name (as a string) from the object. 639*29d8e7bcSMatthew G. Knepley 640*29d8e7bcSMatthew G. Knepley Not Collective 641*29d8e7bcSMatthew G. Knepley 642*29d8e7bcSMatthew G. Knepley Input Parameter: 643*29d8e7bcSMatthew G. Knepley . dm - The PetscSpace 644*29d8e7bcSMatthew G. Knepley 645*29d8e7bcSMatthew G. Knepley Output Parameter: 646*29d8e7bcSMatthew G. Knepley . name - The PetscSpace type name 647*29d8e7bcSMatthew G. Knepley 648*29d8e7bcSMatthew G. Knepley Level: intermediate 649*29d8e7bcSMatthew G. Knepley 650*29d8e7bcSMatthew G. Knepley .keywords: PetscSpace, get, type, name 651*29d8e7bcSMatthew G. Knepley .seealso: PetscSpaceSetType(), PetscSpaceCreate() 652*29d8e7bcSMatthew G. Knepley @*/ 653*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceGetType(PetscSpace sp, PetscSpaceType *name) 654*29d8e7bcSMatthew G. Knepley { 655*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 656*29d8e7bcSMatthew G. Knepley 657*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 658*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 659*29d8e7bcSMatthew G. Knepley PetscValidCharPointer(name, 2); 660*29d8e7bcSMatthew G. Knepley if (!PetscSpaceRegisterAllCalled) { 661*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceRegisterAll();CHKERRQ(ierr); 662*29d8e7bcSMatthew G. Knepley } 663*29d8e7bcSMatthew G. Knepley *name = ((PetscObject) sp)->type_name; 664*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 665*29d8e7bcSMatthew G. Knepley } 666*29d8e7bcSMatthew G. Knepley 667*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 668*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceView" 669*29d8e7bcSMatthew G. Knepley /*@C 670*29d8e7bcSMatthew G. Knepley PetscSpaceView - Views a PetscSpace 671*29d8e7bcSMatthew G. Knepley 672*29d8e7bcSMatthew G. Knepley Collective on PetscSpace 673*29d8e7bcSMatthew G. Knepley 674*29d8e7bcSMatthew G. Knepley Input Parameter: 675*29d8e7bcSMatthew G. Knepley + sp - the PetscSpace object to view 676*29d8e7bcSMatthew G. Knepley - v - the viewer 677*29d8e7bcSMatthew G. Knepley 678*29d8e7bcSMatthew G. Knepley Level: developer 679*29d8e7bcSMatthew G. Knepley 680*29d8e7bcSMatthew G. Knepley .seealso PetscSpaceDestroy() 681*29d8e7bcSMatthew G. Knepley @*/ 682*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceView(PetscSpace sp, PetscViewer v) 683*29d8e7bcSMatthew G. Knepley { 684*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 685*29d8e7bcSMatthew G. Knepley 686*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 687*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 688*29d8e7bcSMatthew G. Knepley if (!v) { 689*29d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr); 690*29d8e7bcSMatthew G. Knepley } 691*29d8e7bcSMatthew G. Knepley if (sp->ops->view) { 692*29d8e7bcSMatthew G. Knepley ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr); 693*29d8e7bcSMatthew G. Knepley } 694*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 695*29d8e7bcSMatthew G. Knepley } 696*29d8e7bcSMatthew G. Knepley 697*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 698*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceViewFromOptions" 699*29d8e7bcSMatthew G. Knepley /* 700*29d8e7bcSMatthew G. Knepley PetscSpaceViewFromOptions - Processes command line options to determine if/how a PetscSpace is to be viewed. 701*29d8e7bcSMatthew G. Knepley 702*29d8e7bcSMatthew G. Knepley Collective on PetscSpace 703*29d8e7bcSMatthew G. Knepley 704*29d8e7bcSMatthew G. Knepley Input Parameters: 705*29d8e7bcSMatthew G. Knepley + sp - the PetscSpace 706*29d8e7bcSMatthew G. Knepley . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd' 707*29d8e7bcSMatthew G. Knepley - optionname - option to activate viewing 708*29d8e7bcSMatthew G. Knepley 709*29d8e7bcSMatthew G. Knepley Level: intermediate 710*29d8e7bcSMatthew G. Knepley 711*29d8e7bcSMatthew G. Knepley .keywords: PetscSpace, view, options, database 712*29d8e7bcSMatthew G. Knepley .seealso: VecViewFromOptions(), MatViewFromOptions() 713*29d8e7bcSMatthew G. Knepley */ 714*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceViewFromOptions(PetscSpace sp, const char prefix[], const char optionname[]) 715*29d8e7bcSMatthew G. Knepley { 716*29d8e7bcSMatthew G. Knepley PetscViewer viewer; 717*29d8e7bcSMatthew G. Knepley PetscViewerFormat format; 718*29d8e7bcSMatthew G. Knepley PetscBool flg; 719*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 720*29d8e7bcSMatthew G. Knepley 721*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 722*29d8e7bcSMatthew G. Knepley if (prefix) { 723*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 724*29d8e7bcSMatthew G. Knepley } else { 725*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 726*29d8e7bcSMatthew G. Knepley } 727*29d8e7bcSMatthew G. Knepley if (flg) { 728*29d8e7bcSMatthew G. Knepley ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 729*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceView(sp, viewer);CHKERRQ(ierr); 730*29d8e7bcSMatthew G. Knepley ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 731*29d8e7bcSMatthew G. Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 732*29d8e7bcSMatthew G. Knepley } 733*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 734*29d8e7bcSMatthew G. Knepley } 735*29d8e7bcSMatthew G. Knepley 736*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 737*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetFromOptions" 738*29d8e7bcSMatthew G. Knepley /*@ 739*29d8e7bcSMatthew G. Knepley PetscSpaceSetFromOptions - sets parameters in a PetscSpace from the options database 740*29d8e7bcSMatthew G. Knepley 741*29d8e7bcSMatthew G. Knepley Collective on PetscSpace 742*29d8e7bcSMatthew G. Knepley 743*29d8e7bcSMatthew G. Knepley Input Parameter: 744*29d8e7bcSMatthew G. Knepley . sp - the PetscSpace object to set options for 745*29d8e7bcSMatthew G. Knepley 746*29d8e7bcSMatthew G. Knepley Options Database: 747*29d8e7bcSMatthew G. Knepley . -petscspace_order the approximation order of the space 748*29d8e7bcSMatthew G. Knepley 749*29d8e7bcSMatthew G. Knepley Level: developer 750*29d8e7bcSMatthew G. Knepley 751*29d8e7bcSMatthew G. Knepley .seealso PetscSpaceView() 752*29d8e7bcSMatthew G. Knepley @*/ 753*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetFromOptions(PetscSpace sp) 754*29d8e7bcSMatthew G. Knepley { 755*29d8e7bcSMatthew G. Knepley const char *defaultType; 756*29d8e7bcSMatthew G. Knepley char typename[256]; 757*29d8e7bcSMatthew G. Knepley PetscBool flg; 758*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 759*29d8e7bcSMatthew G. Knepley 760*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 761*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 762*29d8e7bcSMatthew G. Knepley if (!((PetscObject) sp)->type_name) { 763*29d8e7bcSMatthew G. Knepley defaultType = PETSCSPACEPOLYNOMIAL; 764*29d8e7bcSMatthew G. Knepley } else { 765*29d8e7bcSMatthew G. Knepley defaultType = ((PetscObject) sp)->type_name; 766*29d8e7bcSMatthew G. Knepley } 767*29d8e7bcSMatthew G. Knepley if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 768*29d8e7bcSMatthew G. Knepley 769*29d8e7bcSMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 770*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsList("-petscspace_type", "Linear space", "PetscSpaceSetType", PetscSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr); 771*29d8e7bcSMatthew G. Knepley if (flg) { 772*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceSetType(sp, typename);CHKERRQ(ierr); 773*29d8e7bcSMatthew G. Knepley } else if (!((PetscObject) sp)->type_name) { 774*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceSetType(sp, defaultType);CHKERRQ(ierr); 775*29d8e7bcSMatthew G. Knepley } 776*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsInt("-petscspace_order", "The approximation order", "PetscSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr); 777*29d8e7bcSMatthew G. Knepley if (sp->ops->setfromoptions) { 778*29d8e7bcSMatthew G. Knepley ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr); 779*29d8e7bcSMatthew G. Knepley } 780*29d8e7bcSMatthew G. Knepley /* process any options handlers added with PetscObjectAddOptionsHandler() */ 781*29d8e7bcSMatthew G. Knepley ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr); 782*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 783*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceViewFromOptions(sp, NULL, "-petscspace_view");CHKERRQ(ierr); 784*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 785*29d8e7bcSMatthew G. Knepley } 786*29d8e7bcSMatthew G. Knepley 787*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 788*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetUp" 789*29d8e7bcSMatthew G. Knepley /*@C 790*29d8e7bcSMatthew G. Knepley PetscSpaceSetUp - Construct data structures for the PetscSpace 791*29d8e7bcSMatthew G. Knepley 792*29d8e7bcSMatthew G. Knepley Collective on PetscSpace 793*29d8e7bcSMatthew G. Knepley 794*29d8e7bcSMatthew G. Knepley Input Parameter: 795*29d8e7bcSMatthew G. Knepley . sp - the PetscSpace object to setup 796*29d8e7bcSMatthew G. Knepley 797*29d8e7bcSMatthew G. Knepley Level: developer 798*29d8e7bcSMatthew G. Knepley 799*29d8e7bcSMatthew G. Knepley .seealso PetscSpaceView(), PetscSpaceDestroy() 800*29d8e7bcSMatthew G. Knepley @*/ 801*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetUp(PetscSpace sp) 802*29d8e7bcSMatthew G. Knepley { 803*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 804*29d8e7bcSMatthew G. Knepley 805*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 806*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 807*29d8e7bcSMatthew G. Knepley if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 808*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 809*29d8e7bcSMatthew G. Knepley } 810*29d8e7bcSMatthew G. Knepley 811*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 812*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceDestroy" 813*29d8e7bcSMatthew G. Knepley /*@ 814*29d8e7bcSMatthew G. Knepley PetscSpaceDestroy - Destroys a PetscSpace object 815*29d8e7bcSMatthew G. Knepley 816*29d8e7bcSMatthew G. Knepley Collective on PetscSpace 817*29d8e7bcSMatthew G. Knepley 818*29d8e7bcSMatthew G. Knepley Input Parameter: 819*29d8e7bcSMatthew G. Knepley . sp - the PetscSpace object to destroy 820*29d8e7bcSMatthew G. Knepley 821*29d8e7bcSMatthew G. Knepley Level: developer 822*29d8e7bcSMatthew G. Knepley 823*29d8e7bcSMatthew G. Knepley .seealso PetscSpaceView() 824*29d8e7bcSMatthew G. Knepley @*/ 825*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceDestroy(PetscSpace *sp) 826*29d8e7bcSMatthew G. Knepley { 827*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 828*29d8e7bcSMatthew G. Knepley 829*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 830*29d8e7bcSMatthew G. Knepley if (!*sp) PetscFunctionReturn(0); 831*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific((*sp), PETSCSPACE_CLASSID, 1); 832*29d8e7bcSMatthew G. Knepley 833*29d8e7bcSMatthew G. Knepley if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 834*29d8e7bcSMatthew G. Knepley ((PetscObject) (*sp))->refct = 0; 835*29d8e7bcSMatthew G. Knepley /* if memory was published with AMS then destroy it */ 836*29d8e7bcSMatthew G. Knepley ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr); 837*29d8e7bcSMatthew G. Knepley 838*29d8e7bcSMatthew G. Knepley ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr); 839*29d8e7bcSMatthew G. Knepley 840*29d8e7bcSMatthew G. Knepley ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr); 841*29d8e7bcSMatthew G. Knepley ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 842*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 843*29d8e7bcSMatthew G. Knepley } 844*29d8e7bcSMatthew G. Knepley 845*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 846*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceCreate" 847*29d8e7bcSMatthew G. Knepley /*@ 848*29d8e7bcSMatthew G. Knepley PetscSpaceCreate - Creates an empty PetscSpace object. The type can then be set with PetscSpaceSetType(). 849*29d8e7bcSMatthew G. Knepley 850*29d8e7bcSMatthew G. Knepley Collective on MPI_Comm 851*29d8e7bcSMatthew G. Knepley 852*29d8e7bcSMatthew G. Knepley Input Parameter: 853*29d8e7bcSMatthew G. Knepley . comm - The communicator for the PetscSpace object 854*29d8e7bcSMatthew G. Knepley 855*29d8e7bcSMatthew G. Knepley Output Parameter: 856*29d8e7bcSMatthew G. Knepley . sp - The PetscSpace object 857*29d8e7bcSMatthew G. Knepley 858*29d8e7bcSMatthew G. Knepley Level: beginner 859*29d8e7bcSMatthew G. Knepley 860*29d8e7bcSMatthew G. Knepley .seealso: PetscSpaceSetType(), PETSCSPACEPOLYNOMIAL 861*29d8e7bcSMatthew G. Knepley @*/ 862*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceCreate(MPI_Comm comm, PetscSpace *sp) 863*29d8e7bcSMatthew G. Knepley { 864*29d8e7bcSMatthew G. Knepley PetscSpace s; 865*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 866*29d8e7bcSMatthew G. Knepley 867*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 868*29d8e7bcSMatthew G. Knepley PetscValidPointer(sp, 2); 869*29d8e7bcSMatthew G. Knepley *sp = NULL; 870*29d8e7bcSMatthew G. Knepley #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 871*29d8e7bcSMatthew G. Knepley ierr = PetscFEInitializePackage();CHKERRQ(ierr); 872*29d8e7bcSMatthew G. Knepley #endif 873*29d8e7bcSMatthew G. Knepley 874*29d8e7bcSMatthew G. Knepley ierr = PetscHeaderCreate(s, _p_PetscSpace, struct _PetscSpaceOps, PETSCSPACE_CLASSID, "PetscSpace", "Linear Space", "PetscSpace", comm, PetscSpaceDestroy, PetscSpaceView);CHKERRQ(ierr); 875*29d8e7bcSMatthew G. Knepley ierr = PetscMemzero(s->ops, sizeof(struct _PetscSpaceOps));CHKERRQ(ierr); 876*29d8e7bcSMatthew G. Knepley 877*29d8e7bcSMatthew G. Knepley s->order = 0; 878*29d8e7bcSMatthew G. Knepley ierr = DMShellCreate(comm, &s->dm);CHKERRQ(ierr); 879*29d8e7bcSMatthew G. Knepley 880*29d8e7bcSMatthew G. Knepley *sp = s; 881*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 882*29d8e7bcSMatthew G. Knepley } 883*29d8e7bcSMatthew G. Knepley 884*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 885*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceGetDimension" 886*29d8e7bcSMatthew G. Knepley /* Dimension of the space, i.e. number of basis vectors */ 887*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceGetDimension(PetscSpace sp, PetscInt *dim) 888*29d8e7bcSMatthew G. Knepley { 889*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 890*29d8e7bcSMatthew G. Knepley 891*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 892*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 893*29d8e7bcSMatthew G. Knepley PetscValidPointer(dim, 2); 894*29d8e7bcSMatthew G. Knepley *dim = 0; 895*29d8e7bcSMatthew G. Knepley if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 896*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 897*29d8e7bcSMatthew G. Knepley } 898*29d8e7bcSMatthew G. Knepley 899*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 900*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceGetOrder" 901*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceGetOrder(PetscSpace sp, PetscInt *order) 902*29d8e7bcSMatthew G. Knepley { 903*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 904*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 905*29d8e7bcSMatthew G. Knepley PetscValidPointer(order, 2); 906*29d8e7bcSMatthew G. Knepley *order = sp->order; 907*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 908*29d8e7bcSMatthew G. Knepley } 909*29d8e7bcSMatthew G. Knepley 910*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 911*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetOrder" 912*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetOrder(PetscSpace sp, PetscInt order) 913*29d8e7bcSMatthew G. Knepley { 914*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 915*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 916*29d8e7bcSMatthew G. Knepley sp->order = order; 917*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 918*29d8e7bcSMatthew G. Knepley } 919*29d8e7bcSMatthew G. Knepley 920*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 921*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceEvaluate" 922*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceEvaluate(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 923*29d8e7bcSMatthew G. Knepley { 924*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 925*29d8e7bcSMatthew G. Knepley 926*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 927*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 928*29d8e7bcSMatthew G. Knepley PetscValidPointer(points, 3); 929*29d8e7bcSMatthew G. Knepley if (B) PetscValidPointer(B, 4); 930*29d8e7bcSMatthew G. Knepley if (D) PetscValidPointer(D, 5); 931*29d8e7bcSMatthew G. Knepley if (H) PetscValidPointer(H, 6); 932*29d8e7bcSMatthew G. Knepley if (sp->ops->evaluate) {ierr = (*sp->ops->evaluate)(sp, npoints, points, B, D, H);CHKERRQ(ierr);} 933*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 934*29d8e7bcSMatthew G. Knepley } 935*29d8e7bcSMatthew G. Knepley 936*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 937*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetFromOptions_Polynomial" 938*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscSpace sp) 939*29d8e7bcSMatthew G. Knepley { 940*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 941*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 942*29d8e7bcSMatthew G. Knepley 943*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 944*29d8e7bcSMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 945*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsInt("-petscspace_num_variables", "The number of different variables, e.g. x and y", "PetscSpacePolynomialSetNumVariables", poly->numVariables, &poly->numVariables, NULL);CHKERRQ(ierr); 946*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);CHKERRQ(ierr); 947*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 948*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 949*29d8e7bcSMatthew G. Knepley } 950*29d8e7bcSMatthew G. Knepley 951*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 952*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialView_Ascii" 953*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer viewer) 954*29d8e7bcSMatthew G. Knepley { 955*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 956*29d8e7bcSMatthew G. Knepley PetscViewerFormat format; 957*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 958*29d8e7bcSMatthew G. Knepley 959*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 960*29d8e7bcSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 961*29d8e7bcSMatthew G. Knepley if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 962*29d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr); 963*29d8e7bcSMatthew G. Knepley } else { 964*29d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr); 965*29d8e7bcSMatthew G. Knepley } 966*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 967*29d8e7bcSMatthew G. Knepley } 968*29d8e7bcSMatthew G. Knepley 969*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 970*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceView_Polynomial" 971*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer) 972*29d8e7bcSMatthew G. Knepley { 973*29d8e7bcSMatthew G. Knepley PetscBool iascii; 974*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 975*29d8e7bcSMatthew G. Knepley 976*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 977*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 978*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 979*29d8e7bcSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 980*29d8e7bcSMatthew G. Knepley if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);} 981*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 982*29d8e7bcSMatthew G. Knepley } 983*29d8e7bcSMatthew G. Knepley 984*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 985*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetUp_Polynomial" 986*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp) 987*29d8e7bcSMatthew G. Knepley { 988*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 989*29d8e7bcSMatthew G. Knepley PetscInt ndegree = sp->order+1; 990*29d8e7bcSMatthew G. Knepley PetscInt dim = poly->numVariables; 991*29d8e7bcSMatthew G. Knepley PetscInt pdim, deg; 992*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 993*29d8e7bcSMatthew G. Knepley 994*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 995*29d8e7bcSMatthew G. Knepley ierr = PetscMalloc(ndegree * sizeof(PetscInt), &poly->degrees);CHKERRQ(ierr); 996*29d8e7bcSMatthew G. Knepley for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg; 997*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 998*29d8e7bcSMatthew G. Knepley /* Create A (pdim x ndegree * dim) */ 999*29d8e7bcSMatthew G. Knepley ierr = PetscMalloc(pdim*ndegree*dim * sizeof(PetscReal), &poly->A);CHKERRQ(ierr); 1000*29d8e7bcSMatthew G. Knepley ierr = PetscMemzero(poly->A, pdim*ndegree*dim * sizeof(PetscReal));CHKERRQ(ierr); 1001*29d8e7bcSMatthew G. Knepley /* Hardcode P_1: Here we need a way to iterate through the basis */ 1002*29d8e7bcSMatthew G. Knepley poly->A[(0*ndegree + 0)*dim + 0] = 1.0; /* 1 */ 1003*29d8e7bcSMatthew G. Knepley poly->A[(1*ndegree + 1)*dim + 0] = 1.0; /* x */ 1004*29d8e7bcSMatthew G. Knepley poly->A[(2*ndegree + 1)*dim + 1] = 1.0; /* y */ 1005*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1006*29d8e7bcSMatthew G. Knepley } 1007*29d8e7bcSMatthew G. Knepley 1008*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1009*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceDestroy_Polynomial" 1010*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp) 1011*29d8e7bcSMatthew G. Knepley { 1012*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1013*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1014*29d8e7bcSMatthew G. Knepley 1015*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1016*29d8e7bcSMatthew G. Knepley ierr = PetscFree(poly->A);CHKERRQ(ierr); 1017*29d8e7bcSMatthew G. Knepley ierr = PetscFree(poly->degrees);CHKERRQ(ierr); 1018*29d8e7bcSMatthew G. Knepley ierr = PetscFree(poly);CHKERRQ(ierr); 1019*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1020*29d8e7bcSMatthew G. Knepley } 1021*29d8e7bcSMatthew G. Knepley 1022*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1023*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceGetDimension_Polynomial" 1024*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim) 1025*29d8e7bcSMatthew G. Knepley { 1026*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1027*29d8e7bcSMatthew G. Knepley PetscInt deg = sp->order; 1028*29d8e7bcSMatthew G. Knepley PetscInt n = poly->numVariables, i; 1029*29d8e7bcSMatthew G. Knepley PetscReal D = 1.0; 1030*29d8e7bcSMatthew G. Knepley 1031*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1032*29d8e7bcSMatthew G. Knepley for (i = 1; i <= n; ++i) { 1033*29d8e7bcSMatthew G. Knepley D *= ((PetscReal) (deg+i))/i; 1034*29d8e7bcSMatthew G. Knepley } 1035*29d8e7bcSMatthew G. Knepley *dim = (PetscInt) (D + 0.5); 1036*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1037*29d8e7bcSMatthew G. Knepley } 1038*29d8e7bcSMatthew G. Knepley 1039*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1040*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceEvaluate_Polynomial" 1041*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 1042*29d8e7bcSMatthew G. Knepley { 1043*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1044*29d8e7bcSMatthew G. Knepley DM dm = sp->dm; 1045*29d8e7bcSMatthew G. Knepley PetscInt ndegree = sp->order+1; 1046*29d8e7bcSMatthew G. Knepley PetscInt *degrees = poly->degrees; 1047*29d8e7bcSMatthew G. Knepley PetscInt dim = poly->numVariables; 1048*29d8e7bcSMatthew G. Knepley PetscReal *A = poly->A; 1049*29d8e7bcSMatthew G. Knepley PetscReal *lpoints, *tmp, *LB, *LD, *LH; 1050*29d8e7bcSMatthew G. Knepley PetscInt pdim, d, i, p, deg; 1051*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1052*29d8e7bcSMatthew G. Knepley 1053*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1054*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 1055*29d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr); 1056*29d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr); 1057*29d8e7bcSMatthew G. Knepley if (B) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);} 1058*29d8e7bcSMatthew G. Knepley if (D) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);} 1059*29d8e7bcSMatthew G. Knepley if (H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);} 1060*29d8e7bcSMatthew G. Knepley for (d = 0; d < dim; ++d) { 1061*29d8e7bcSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 1062*29d8e7bcSMatthew G. Knepley lpoints[p] = points[p*dim+d]; 1063*29d8e7bcSMatthew G. Knepley } 1064*29d8e7bcSMatthew G. Knepley ierr = PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);CHKERRQ(ierr); 1065*29d8e7bcSMatthew G. Knepley /* LB, LD, LH (ndegree * dim x npoints) */ 1066*29d8e7bcSMatthew G. Knepley for (deg = 0; deg < ndegree; ++deg) { 1067*29d8e7bcSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 1068*29d8e7bcSMatthew G. Knepley if (B) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg]; 1069*29d8e7bcSMatthew G. Knepley if (D) LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg]; 1070*29d8e7bcSMatthew G. Knepley if (H) LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg]; 1071*29d8e7bcSMatthew G. Knepley } 1072*29d8e7bcSMatthew G. Knepley } 1073*29d8e7bcSMatthew G. Knepley } 1074*29d8e7bcSMatthew G. Knepley /* Multiply by A (pdim x ndegree * dim) */ 1075*29d8e7bcSMatthew G. Knepley if (B) { 1076*29d8e7bcSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 1077*29d8e7bcSMatthew G. Knepley for (i = 0; i < pdim; ++i) { 1078*29d8e7bcSMatthew G. Knepley B[p*pdim + i] = 0.0; 1079*29d8e7bcSMatthew G. Knepley for (deg = 0; deg < ndegree; ++deg) { 1080*29d8e7bcSMatthew G. Knepley for (d = 0; d < dim; ++d) { 1081*29d8e7bcSMatthew G. Knepley B[p*pdim + i] += A[(i*ndegree + deg)*dim + d] * LB[(deg*dim + d)*npoints + p]; 1082*29d8e7bcSMatthew G. Knepley } 1083*29d8e7bcSMatthew G. Knepley } 1084*29d8e7bcSMatthew G. Knepley } 1085*29d8e7bcSMatthew G. Knepley } 1086*29d8e7bcSMatthew G. Knepley } 1087*29d8e7bcSMatthew G. Knepley if (D) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code first derivatives"); 1088*29d8e7bcSMatthew G. Knepley if (H) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code second derivatives"); 1089*29d8e7bcSMatthew G. Knepley if (B) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);} 1090*29d8e7bcSMatthew G. Knepley if (D) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);} 1091*29d8e7bcSMatthew G. Knepley if (H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);} 1092*29d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr); 1093*29d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr); 1094*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1095*29d8e7bcSMatthew G. Knepley } 1096*29d8e7bcSMatthew G. Knepley 1097*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1098*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceInitialize_Polynomial" 1099*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp) 1100*29d8e7bcSMatthew G. Knepley { 1101*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1102*29d8e7bcSMatthew G. Knepley sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial; 1103*29d8e7bcSMatthew G. Knepley sp->ops->setup = PetscSpaceSetUp_Polynomial; 1104*29d8e7bcSMatthew G. Knepley sp->ops->view = PetscSpaceView_Polynomial; 1105*29d8e7bcSMatthew G. Knepley sp->ops->destroy = PetscSpaceDestroy_Polynomial; 1106*29d8e7bcSMatthew G. Knepley sp->ops->getdimension = PetscSpaceGetDimension_Polynomial; 1107*29d8e7bcSMatthew G. Knepley sp->ops->evaluate = PetscSpaceEvaluate_Polynomial; 1108*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1109*29d8e7bcSMatthew G. Knepley } 1110*29d8e7bcSMatthew G. Knepley 1111*29d8e7bcSMatthew G. Knepley /*MC 1112*29d8e7bcSMatthew G. Knepley PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of linear polynomials. 1113*29d8e7bcSMatthew G. Knepley 1114*29d8e7bcSMatthew G. Knepley Level: intermediate 1115*29d8e7bcSMatthew G. Knepley 1116*29d8e7bcSMatthew G. Knepley .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType() 1117*29d8e7bcSMatthew G. Knepley M*/ 1118*29d8e7bcSMatthew G. Knepley 1119*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1120*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceCreate_Polynomial" 1121*29d8e7bcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp) 1122*29d8e7bcSMatthew G. Knepley { 1123*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly; 1124*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1125*29d8e7bcSMatthew G. Knepley 1126*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1127*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1128*29d8e7bcSMatthew G. Knepley ierr = PetscNewLog(sp, PetscSpace_Poly, &poly);CHKERRQ(ierr); 1129*29d8e7bcSMatthew G. Knepley sp->data = poly; 1130*29d8e7bcSMatthew G. Knepley 1131*29d8e7bcSMatthew G. Knepley poly->numVariables = 0; 1132*29d8e7bcSMatthew G. Knepley poly->symmetric = PETSC_FALSE; 1133*29d8e7bcSMatthew G. Knepley poly->degrees = NULL; 1134*29d8e7bcSMatthew G. Knepley poly->A = NULL; 1135*29d8e7bcSMatthew G. Knepley 1136*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr); 1137*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1138*29d8e7bcSMatthew G. Knepley } 1139*29d8e7bcSMatthew G. Knepley 1140*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1141*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialSetSymmetric" 1142*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool sym) 1143*29d8e7bcSMatthew G. Knepley { 1144*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1145*29d8e7bcSMatthew G. Knepley 1146*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1147*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1148*29d8e7bcSMatthew G. Knepley poly->symmetric = sym; 1149*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1150*29d8e7bcSMatthew G. Knepley } 1151*29d8e7bcSMatthew G. Knepley 1152*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1153*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialGetSymmetric" 1154*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace sp, PetscBool *sym) 1155*29d8e7bcSMatthew G. Knepley { 1156*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1157*29d8e7bcSMatthew G. Knepley 1158*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1159*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1160*29d8e7bcSMatthew G. Knepley PetscValidPointer(sym, 2); 1161*29d8e7bcSMatthew G. Knepley *sym = poly->symmetric; 1162*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1163*29d8e7bcSMatthew G. Knepley } 1164*29d8e7bcSMatthew G. Knepley 1165*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1166*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialSetNumVariables" 1167*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace sp, PetscInt n) 1168*29d8e7bcSMatthew G. Knepley { 1169*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1170*29d8e7bcSMatthew G. Knepley 1171*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1172*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1173*29d8e7bcSMatthew G. Knepley poly->numVariables = n; 1174*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1175*29d8e7bcSMatthew G. Knepley } 1176*29d8e7bcSMatthew G. Knepley 1177*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1178*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialGetNumVariables" 1179*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace sp, PetscInt *n) 1180*29d8e7bcSMatthew G. Knepley { 1181*29d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 1182*29d8e7bcSMatthew G. Knepley 1183*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1184*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 1185*29d8e7bcSMatthew G. Knepley PetscValidPointer(n, 2); 1186*29d8e7bcSMatthew G. Knepley *n = poly->numVariables; 1187*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1188*29d8e7bcSMatthew G. Knepley } 1189*29d8e7bcSMatthew G. Knepley 1190*29d8e7bcSMatthew G. Knepley 1191*29d8e7bcSMatthew G. Knepley PetscInt PETSCDUALSPACE_CLASSID = 0; 1192*29d8e7bcSMatthew G. Knepley 1193*29d8e7bcSMatthew G. Knepley PetscFunctionList PetscDualSpaceList = NULL; 1194*29d8e7bcSMatthew G. Knepley PetscBool PetscDualSpaceRegisterAllCalled = PETSC_FALSE; 1195*29d8e7bcSMatthew G. Knepley 1196*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1197*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceRegister" 1198*29d8e7bcSMatthew G. Knepley /*@C 1199*29d8e7bcSMatthew G. Knepley PetscDualSpaceRegister - Adds a new PetscDualSpace implementation 1200*29d8e7bcSMatthew G. Knepley 1201*29d8e7bcSMatthew G. Knepley Not Collective 1202*29d8e7bcSMatthew G. Knepley 1203*29d8e7bcSMatthew G. Knepley Input Parameters: 1204*29d8e7bcSMatthew G. Knepley + name - The name of a new user-defined creation routine 1205*29d8e7bcSMatthew G. Knepley - create_func - The creation routine itself 1206*29d8e7bcSMatthew G. Knepley 1207*29d8e7bcSMatthew G. Knepley Notes: 1208*29d8e7bcSMatthew G. Knepley PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces 1209*29d8e7bcSMatthew G. Knepley 1210*29d8e7bcSMatthew G. Knepley Sample usage: 1211*29d8e7bcSMatthew G. Knepley .vb 1212*29d8e7bcSMatthew G. Knepley PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate); 1213*29d8e7bcSMatthew G. Knepley .ve 1214*29d8e7bcSMatthew G. Knepley 1215*29d8e7bcSMatthew G. Knepley Then, your PetscDualSpace type can be chosen with the procedural interface via 1216*29d8e7bcSMatthew G. Knepley .vb 1217*29d8e7bcSMatthew G. Knepley PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 1218*29d8e7bcSMatthew G. Knepley PetscDualSpaceSetType(PetscDualSpace, "my_dual_space"); 1219*29d8e7bcSMatthew G. Knepley .ve 1220*29d8e7bcSMatthew G. Knepley or at runtime via the option 1221*29d8e7bcSMatthew G. Knepley .vb 1222*29d8e7bcSMatthew G. Knepley -petscdualspace_type my_dual_space 1223*29d8e7bcSMatthew G. Knepley .ve 1224*29d8e7bcSMatthew G. Knepley 1225*29d8e7bcSMatthew G. Knepley Level: advanced 1226*29d8e7bcSMatthew G. Knepley 1227*29d8e7bcSMatthew G. Knepley .keywords: PetscDualSpace, register 1228*29d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy() 1229*29d8e7bcSMatthew G. Knepley 1230*29d8e7bcSMatthew G. Knepley @*/ 1231*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace)) 1232*29d8e7bcSMatthew G. Knepley { 1233*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1234*29d8e7bcSMatthew G. Knepley 1235*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1236*29d8e7bcSMatthew G. Knepley ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr); 1237*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1238*29d8e7bcSMatthew G. Knepley } 1239*29d8e7bcSMatthew G. Knepley 1240*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1241*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetType" 1242*29d8e7bcSMatthew G. Knepley /*@C 1243*29d8e7bcSMatthew G. Knepley PetscDualSpaceSetType - Builds a particular PetscDualSpace 1244*29d8e7bcSMatthew G. Knepley 1245*29d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 1246*29d8e7bcSMatthew G. Knepley 1247*29d8e7bcSMatthew G. Knepley Input Parameters: 1248*29d8e7bcSMatthew G. Knepley + sp - The PetscDualSpace object 1249*29d8e7bcSMatthew G. Knepley - name - The kind of space 1250*29d8e7bcSMatthew G. Knepley 1251*29d8e7bcSMatthew G. Knepley Options Database Key: 1252*29d8e7bcSMatthew G. Knepley . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types 1253*29d8e7bcSMatthew G. Knepley 1254*29d8e7bcSMatthew G. Knepley Level: intermediate 1255*29d8e7bcSMatthew G. Knepley 1256*29d8e7bcSMatthew G. Knepley .keywords: PetscDualSpace, set, type 1257*29d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate() 1258*29d8e7bcSMatthew G. Knepley @*/ 1259*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name) 1260*29d8e7bcSMatthew G. Knepley { 1261*29d8e7bcSMatthew G. Knepley PetscErrorCode (*r)(PetscDualSpace); 1262*29d8e7bcSMatthew G. Knepley PetscBool match; 1263*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1264*29d8e7bcSMatthew G. Knepley 1265*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1266*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1267*29d8e7bcSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 1268*29d8e7bcSMatthew G. Knepley if (match) PetscFunctionReturn(0); 1269*29d8e7bcSMatthew G. Knepley 1270*29d8e7bcSMatthew G. Knepley if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);} 1271*29d8e7bcSMatthew G. Knepley ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr); 1272*29d8e7bcSMatthew G. Knepley if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name); 1273*29d8e7bcSMatthew G. Knepley 1274*29d8e7bcSMatthew G. Knepley if (sp->ops->destroy) { 1275*29d8e7bcSMatthew G. Knepley ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 1276*29d8e7bcSMatthew G. Knepley sp->ops->destroy = NULL; 1277*29d8e7bcSMatthew G. Knepley } 1278*29d8e7bcSMatthew G. Knepley ierr = (*r)(sp);CHKERRQ(ierr); 1279*29d8e7bcSMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 1280*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1281*29d8e7bcSMatthew G. Knepley } 1282*29d8e7bcSMatthew G. Knepley 1283*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1284*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetType" 1285*29d8e7bcSMatthew G. Knepley /*@C 1286*29d8e7bcSMatthew G. Knepley PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object. 1287*29d8e7bcSMatthew G. Knepley 1288*29d8e7bcSMatthew G. Knepley Not Collective 1289*29d8e7bcSMatthew G. Knepley 1290*29d8e7bcSMatthew G. Knepley Input Parameter: 1291*29d8e7bcSMatthew G. Knepley . dm - The PetscDualSpace 1292*29d8e7bcSMatthew G. Knepley 1293*29d8e7bcSMatthew G. Knepley Output Parameter: 1294*29d8e7bcSMatthew G. Knepley . name - The PetscDualSpace type name 1295*29d8e7bcSMatthew G. Knepley 1296*29d8e7bcSMatthew G. Knepley Level: intermediate 1297*29d8e7bcSMatthew G. Knepley 1298*29d8e7bcSMatthew G. Knepley .keywords: PetscDualSpace, get, type, name 1299*29d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate() 1300*29d8e7bcSMatthew G. Knepley @*/ 1301*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name) 1302*29d8e7bcSMatthew G. Knepley { 1303*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1304*29d8e7bcSMatthew G. Knepley 1305*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1306*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1307*29d8e7bcSMatthew G. Knepley PetscValidCharPointer(name, 2); 1308*29d8e7bcSMatthew G. Knepley if (!PetscDualSpaceRegisterAllCalled) { 1309*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr); 1310*29d8e7bcSMatthew G. Knepley } 1311*29d8e7bcSMatthew G. Knepley *name = ((PetscObject) sp)->type_name; 1312*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1313*29d8e7bcSMatthew G. Knepley } 1314*29d8e7bcSMatthew G. Knepley 1315*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1316*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceView" 1317*29d8e7bcSMatthew G. Knepley /*@C 1318*29d8e7bcSMatthew G. Knepley PetscDualSpaceView - Views a PetscDualSpace 1319*29d8e7bcSMatthew G. Knepley 1320*29d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 1321*29d8e7bcSMatthew G. Knepley 1322*29d8e7bcSMatthew G. Knepley Input Parameter: 1323*29d8e7bcSMatthew G. Knepley + sp - the PetscDualSpace object to view 1324*29d8e7bcSMatthew G. Knepley - v - the viewer 1325*29d8e7bcSMatthew G. Knepley 1326*29d8e7bcSMatthew G. Knepley Level: developer 1327*29d8e7bcSMatthew G. Knepley 1328*29d8e7bcSMatthew G. Knepley .seealso PetscDualSpaceDestroy() 1329*29d8e7bcSMatthew G. Knepley @*/ 1330*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v) 1331*29d8e7bcSMatthew G. Knepley { 1332*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1333*29d8e7bcSMatthew G. Knepley 1334*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1335*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1336*29d8e7bcSMatthew G. Knepley if (!v) { 1337*29d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr); 1338*29d8e7bcSMatthew G. Knepley } 1339*29d8e7bcSMatthew G. Knepley if (sp->ops->view) { 1340*29d8e7bcSMatthew G. Knepley ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr); 1341*29d8e7bcSMatthew G. Knepley } 1342*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1343*29d8e7bcSMatthew G. Knepley } 1344*29d8e7bcSMatthew G. Knepley 1345*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1346*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceViewFromOptions" 1347*29d8e7bcSMatthew G. Knepley /* 1348*29d8e7bcSMatthew G. Knepley PetscDualSpaceViewFromOptions - Processes command line options to determine if/how a PetscDualSpace is to be viewed. 1349*29d8e7bcSMatthew G. Knepley 1350*29d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 1351*29d8e7bcSMatthew G. Knepley 1352*29d8e7bcSMatthew G. Knepley Input Parameters: 1353*29d8e7bcSMatthew G. Knepley + sp - the PetscDualSpace 1354*29d8e7bcSMatthew G. Knepley . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd' 1355*29d8e7bcSMatthew G. Knepley - optionname - option to activate viewing 1356*29d8e7bcSMatthew G. Knepley 1357*29d8e7bcSMatthew G. Knepley Level: intermediate 1358*29d8e7bcSMatthew G. Knepley 1359*29d8e7bcSMatthew G. Knepley .keywords: PetscDualSpace, view, options, database 1360*29d8e7bcSMatthew G. Knepley .seealso: VecViewFromOptions(), MatViewFromOptions() 1361*29d8e7bcSMatthew G. Knepley */ 1362*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace sp, const char prefix[], const char optionname[]) 1363*29d8e7bcSMatthew G. Knepley { 1364*29d8e7bcSMatthew G. Knepley PetscViewer viewer; 1365*29d8e7bcSMatthew G. Knepley PetscViewerFormat format; 1366*29d8e7bcSMatthew G. Knepley PetscBool flg; 1367*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1368*29d8e7bcSMatthew G. Knepley 1369*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1370*29d8e7bcSMatthew G. Knepley if (prefix) { 1371*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 1372*29d8e7bcSMatthew G. Knepley } else { 1373*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 1374*29d8e7bcSMatthew G. Knepley } 1375*29d8e7bcSMatthew G. Knepley if (flg) { 1376*29d8e7bcSMatthew G. Knepley ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 1377*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceView(sp, viewer);CHKERRQ(ierr); 1378*29d8e7bcSMatthew G. Knepley ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1379*29d8e7bcSMatthew G. Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1380*29d8e7bcSMatthew G. Knepley } 1381*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1382*29d8e7bcSMatthew G. Knepley } 1383*29d8e7bcSMatthew G. Knepley 1384*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1385*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetFromOptions" 1386*29d8e7bcSMatthew G. Knepley /*@ 1387*29d8e7bcSMatthew G. Knepley PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database 1388*29d8e7bcSMatthew G. Knepley 1389*29d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 1390*29d8e7bcSMatthew G. Knepley 1391*29d8e7bcSMatthew G. Knepley Input Parameter: 1392*29d8e7bcSMatthew G. Knepley . sp - the PetscDualSpace object to set options for 1393*29d8e7bcSMatthew G. Knepley 1394*29d8e7bcSMatthew G. Knepley Options Database: 1395*29d8e7bcSMatthew G. Knepley . -petscspace_order the approximation order of the space 1396*29d8e7bcSMatthew G. Knepley 1397*29d8e7bcSMatthew G. Knepley Level: developer 1398*29d8e7bcSMatthew G. Knepley 1399*29d8e7bcSMatthew G. Knepley .seealso PetscDualSpaceView() 1400*29d8e7bcSMatthew G. Knepley @*/ 1401*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp) 1402*29d8e7bcSMatthew G. Knepley { 1403*29d8e7bcSMatthew G. Knepley const char *defaultType; 1404*29d8e7bcSMatthew G. Knepley char typename[256]; 1405*29d8e7bcSMatthew G. Knepley PetscBool flg; 1406*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1407*29d8e7bcSMatthew G. Knepley 1408*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1409*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1410*29d8e7bcSMatthew G. Knepley if (!((PetscObject) sp)->type_name) { 1411*29d8e7bcSMatthew G. Knepley defaultType = PETSCDUALSPACELAGRANGE; 1412*29d8e7bcSMatthew G. Knepley } else { 1413*29d8e7bcSMatthew G. Knepley defaultType = ((PetscObject) sp)->type_name; 1414*29d8e7bcSMatthew G. Knepley } 1415*29d8e7bcSMatthew G. Knepley if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 1416*29d8e7bcSMatthew G. Knepley 1417*29d8e7bcSMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 1418*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr); 1419*29d8e7bcSMatthew G. Knepley if (flg) { 1420*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceSetType(sp, typename);CHKERRQ(ierr); 1421*29d8e7bcSMatthew G. Knepley } else if (!((PetscObject) sp)->type_name) { 1422*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr); 1423*29d8e7bcSMatthew G. Knepley } 1424*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr); 1425*29d8e7bcSMatthew G. Knepley if (sp->ops->setfromoptions) { 1426*29d8e7bcSMatthew G. Knepley ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr); 1427*29d8e7bcSMatthew G. Knepley } 1428*29d8e7bcSMatthew G. Knepley /* process any options handlers added with PetscObjectAddOptionsHandler() */ 1429*29d8e7bcSMatthew G. Knepley ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr); 1430*29d8e7bcSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 1431*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr); 1432*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1433*29d8e7bcSMatthew G. Knepley } 1434*29d8e7bcSMatthew G. Knepley 1435*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1436*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetUp" 1437*29d8e7bcSMatthew G. Knepley /*@C 1438*29d8e7bcSMatthew G. Knepley PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace 1439*29d8e7bcSMatthew G. Knepley 1440*29d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 1441*29d8e7bcSMatthew G. Knepley 1442*29d8e7bcSMatthew G. Knepley Input Parameter: 1443*29d8e7bcSMatthew G. Knepley . sp - the PetscDualSpace object to setup 1444*29d8e7bcSMatthew G. Knepley 1445*29d8e7bcSMatthew G. Knepley Level: developer 1446*29d8e7bcSMatthew G. Knepley 1447*29d8e7bcSMatthew G. Knepley .seealso PetscDualSpaceView(), PetscDualSpaceDestroy() 1448*29d8e7bcSMatthew G. Knepley @*/ 1449*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp) 1450*29d8e7bcSMatthew G. Knepley { 1451*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1452*29d8e7bcSMatthew G. Knepley 1453*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1454*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1455*29d8e7bcSMatthew G. Knepley if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 1456*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1457*29d8e7bcSMatthew G. Knepley } 1458*29d8e7bcSMatthew G. Knepley 1459*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1460*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceDestroy" 1461*29d8e7bcSMatthew G. Knepley /*@ 1462*29d8e7bcSMatthew G. Knepley PetscDualSpaceDestroy - Destroys a PetscDualSpace object 1463*29d8e7bcSMatthew G. Knepley 1464*29d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 1465*29d8e7bcSMatthew G. Knepley 1466*29d8e7bcSMatthew G. Knepley Input Parameter: 1467*29d8e7bcSMatthew G. Knepley . sp - the PetscDualSpace object to destroy 1468*29d8e7bcSMatthew G. Knepley 1469*29d8e7bcSMatthew G. Knepley Level: developer 1470*29d8e7bcSMatthew G. Knepley 1471*29d8e7bcSMatthew G. Knepley .seealso PetscDualSpaceView() 1472*29d8e7bcSMatthew G. Knepley @*/ 1473*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp) 1474*29d8e7bcSMatthew G. Knepley { 1475*29d8e7bcSMatthew G. Knepley PetscInt dim, f; 1476*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1477*29d8e7bcSMatthew G. Knepley 1478*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1479*29d8e7bcSMatthew G. Knepley if (!*sp) PetscFunctionReturn(0); 1480*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1); 1481*29d8e7bcSMatthew G. Knepley 1482*29d8e7bcSMatthew G. Knepley if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 1483*29d8e7bcSMatthew G. Knepley ((PetscObject) (*sp))->refct = 0; 1484*29d8e7bcSMatthew G. Knepley /* if memory was published with AMS then destroy it */ 1485*29d8e7bcSMatthew G. Knepley ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr); 1486*29d8e7bcSMatthew G. Knepley 1487*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr); 1488*29d8e7bcSMatthew G. Knepley for (f = 0; f < dim; ++f) { 1489*29d8e7bcSMatthew G. Knepley /* ierr = PetscQuadratureDestroy((*sp)->functional[f]);CHKERRQ(ierr); */ 1490*29d8e7bcSMatthew G. Knepley } 1491*29d8e7bcSMatthew G. Knepley ierr = PetscFree((*sp)->functional);CHKERRQ(ierr); 1492*29d8e7bcSMatthew G. Knepley 1493*29d8e7bcSMatthew G. Knepley ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr); 1494*29d8e7bcSMatthew G. Knepley ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 1495*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1496*29d8e7bcSMatthew G. Knepley } 1497*29d8e7bcSMatthew G. Knepley 1498*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1499*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceCreate" 1500*29d8e7bcSMatthew G. Knepley /*@ 1501*29d8e7bcSMatthew G. Knepley PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType(). 1502*29d8e7bcSMatthew G. Knepley 1503*29d8e7bcSMatthew G. Knepley Collective on MPI_Comm 1504*29d8e7bcSMatthew G. Knepley 1505*29d8e7bcSMatthew G. Knepley Input Parameter: 1506*29d8e7bcSMatthew G. Knepley . comm - The communicator for the PetscDualSpace object 1507*29d8e7bcSMatthew G. Knepley 1508*29d8e7bcSMatthew G. Knepley Output Parameter: 1509*29d8e7bcSMatthew G. Knepley . sp - The PetscDualSpace object 1510*29d8e7bcSMatthew G. Knepley 1511*29d8e7bcSMatthew G. Knepley Level: beginner 1512*29d8e7bcSMatthew G. Knepley 1513*29d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE 1514*29d8e7bcSMatthew G. Knepley @*/ 1515*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp) 1516*29d8e7bcSMatthew G. Knepley { 1517*29d8e7bcSMatthew G. Knepley PetscDualSpace s; 1518*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1519*29d8e7bcSMatthew G. Knepley 1520*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1521*29d8e7bcSMatthew G. Knepley PetscValidPointer(sp, 2); 1522*29d8e7bcSMatthew G. Knepley *sp = NULL; 1523*29d8e7bcSMatthew G. Knepley #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1524*29d8e7bcSMatthew G. Knepley ierr = PetscFEInitializePackage();CHKERRQ(ierr); 1525*29d8e7bcSMatthew G. Knepley #endif 1526*29d8e7bcSMatthew G. Knepley 1527*29d8e7bcSMatthew G. Knepley ierr = PetscHeaderCreate(s, _p_PetscDualSpace, struct _PetscDualSpaceOps, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr); 1528*29d8e7bcSMatthew G. Knepley ierr = PetscMemzero(s->ops, sizeof(struct _PetscDualSpaceOps));CHKERRQ(ierr); 1529*29d8e7bcSMatthew G. Knepley 1530*29d8e7bcSMatthew G. Knepley s->order = 0; 1531*29d8e7bcSMatthew G. Knepley 1532*29d8e7bcSMatthew G. Knepley *sp = s; 1533*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1534*29d8e7bcSMatthew G. Knepley } 1535*29d8e7bcSMatthew G. Knepley 1536*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1537*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetDM" 1538*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm) 1539*29d8e7bcSMatthew G. Knepley { 1540*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1541*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1542*29d8e7bcSMatthew G. Knepley PetscValidPointer(dm, 2); 1543*29d8e7bcSMatthew G. Knepley *dm = sp->dm; 1544*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1545*29d8e7bcSMatthew G. Knepley } 1546*29d8e7bcSMatthew G. Knepley 1547*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1548*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetDM" 1549*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm) 1550*29d8e7bcSMatthew G. Knepley { 1551*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1552*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1553*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1554*29d8e7bcSMatthew G. Knepley sp->dm = dm; 1555*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1556*29d8e7bcSMatthew G. Knepley } 1557*29d8e7bcSMatthew G. Knepley 1558*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1559*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetOrder" 1560*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order) 1561*29d8e7bcSMatthew G. Knepley { 1562*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1563*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1564*29d8e7bcSMatthew G. Knepley PetscValidPointer(order, 2); 1565*29d8e7bcSMatthew G. Knepley *order = sp->order; 1566*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1567*29d8e7bcSMatthew G. Knepley } 1568*29d8e7bcSMatthew G. Knepley 1569*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1570*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetOrder" 1571*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order) 1572*29d8e7bcSMatthew G. Knepley { 1573*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1574*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1575*29d8e7bcSMatthew G. Knepley sp->order = order; 1576*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1577*29d8e7bcSMatthew G. Knepley } 1578*29d8e7bcSMatthew G. Knepley 1579*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1580*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetFunctional" 1581*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional) 1582*29d8e7bcSMatthew G. Knepley { 1583*29d8e7bcSMatthew G. Knepley PetscInt dim; 1584*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1585*29d8e7bcSMatthew G. Knepley 1586*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1587*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1588*29d8e7bcSMatthew G. Knepley PetscValidPointer(functional, 3); 1589*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 1590*29d8e7bcSMatthew G. Knepley if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim); 1591*29d8e7bcSMatthew G. Knepley *functional = sp->functional[i]; 1592*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1593*29d8e7bcSMatthew G. Knepley } 1594*29d8e7bcSMatthew G. Knepley 1595*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1596*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetDimension" 1597*29d8e7bcSMatthew G. Knepley /* Dimension of the space, i.e. number of basis vectors */ 1598*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim) 1599*29d8e7bcSMatthew G. Knepley { 1600*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1601*29d8e7bcSMatthew G. Knepley 1602*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1603*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1604*29d8e7bcSMatthew G. Knepley PetscValidPointer(dim, 2); 1605*29d8e7bcSMatthew G. Knepley *dim = 0; 1606*29d8e7bcSMatthew G. Knepley if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 1607*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1608*29d8e7bcSMatthew G. Knepley } 1609*29d8e7bcSMatthew G. Knepley 1610*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1611*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetUp_Lagrange" 1612*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 1613*29d8e7bcSMatthew G. Knepley { 1614*29d8e7bcSMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 1615*29d8e7bcSMatthew G. Knepley DM dm = sp->dm; 1616*29d8e7bcSMatthew G. Knepley PetscInt order = sp->order; 1617*29d8e7bcSMatthew G. Knepley PetscSection csection; 1618*29d8e7bcSMatthew G. Knepley Vec coordinates; 1619*29d8e7bcSMatthew G. Knepley PetscScalar *coords; 1620*29d8e7bcSMatthew G. Knepley PetscInt pdim, dim, vStart, vEnd, cStart, coneSize, f = 0; 1621*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1622*29d8e7bcSMatthew G. Knepley 1623*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1624*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 1625*29d8e7bcSMatthew G. Knepley ierr = PetscMalloc(pdim * sizeof(PetscQuadrature), &sp->functional);CHKERRQ(ierr); 1626*29d8e7bcSMatthew G. Knepley /* Classify element type */ 1627*29d8e7bcSMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1628*29d8e7bcSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 1629*29d8e7bcSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1630*29d8e7bcSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cStart, &coneSize);CHKERRQ(ierr); 1631*29d8e7bcSMatthew G. Knepley ierr = DMPlexGetCoordinateSection(dm, &csection);CHKERRQ(ierr); 1632*29d8e7bcSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1633*29d8e7bcSMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1634*29d8e7bcSMatthew G. Knepley if (coneSize == dim+1) { 1635*29d8e7bcSMatthew G. Knepley PetscInt *closure = NULL, closureSize, c; 1636*29d8e7bcSMatthew G. Knepley 1637*29d8e7bcSMatthew G. Knepley /* Simplex */ 1638*29d8e7bcSMatthew G. Knepley if (order > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle order %d basis", order); 1639*29d8e7bcSMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cStart, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1640*29d8e7bcSMatthew G. Knepley for (c = 0; c < closureSize*2; c += 2) { 1641*29d8e7bcSMatthew G. Knepley const PetscInt p = closure[c]; 1642*29d8e7bcSMatthew G. Knepley /* Vertices */ 1643*29d8e7bcSMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) { 1644*29d8e7bcSMatthew G. Knepley PetscReal *qpoints, *qweights; 1645*29d8e7bcSMatthew G. Knepley PetscInt dof, off, d; 1646*29d8e7bcSMatthew G. Knepley 1647*29d8e7bcSMatthew G. Knepley sp->functional[f].numQuadPoints = 1; 1648*29d8e7bcSMatthew G. Knepley ierr = PetscMalloc(sp->functional[f].numQuadPoints*dim * sizeof(PetscReal), &qpoints);CHKERRQ(ierr); 1649*29d8e7bcSMatthew G. Knepley ierr = PetscMalloc(sp->functional[f].numQuadPoints * sizeof(PetscReal), &qweights);CHKERRQ(ierr); 1650*29d8e7bcSMatthew G. Knepley ierr = PetscSectionGetDof(csection, p, &dof);CHKERRQ(ierr); 1651*29d8e7bcSMatthew G. Knepley ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr); 1652*29d8e7bcSMatthew G. Knepley if (dof != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of coordinates %d does not match spatial dimension %d", dof, dim); 1653*29d8e7bcSMatthew G. Knepley for (d = 0; d < dof; ++d) {qpoints[d] = coords[off+d];} 1654*29d8e7bcSMatthew G. Knepley qweights[0] = 1.0; 1655*29d8e7bcSMatthew G. Knepley sp->functional[f].quadPoints = qpoints; 1656*29d8e7bcSMatthew G. Knepley sp->functional[f].quadWeights = qweights; 1657*29d8e7bcSMatthew G. Knepley ++f; 1658*29d8e7bcSMatthew G. Knepley } 1659*29d8e7bcSMatthew G. Knepley } 1660*29d8e7bcSMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cStart, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1661*29d8e7bcSMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle cells with cone size %d", coneSize); 1662*29d8e7bcSMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1663*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1664*29d8e7bcSMatthew G. Knepley } 1665*29d8e7bcSMatthew G. Knepley 1666*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1667*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetDimension_Lagrange" 1668*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim) 1669*29d8e7bcSMatthew G. Knepley { 1670*29d8e7bcSMatthew G. Knepley PetscInt deg = sp->order; 1671*29d8e7bcSMatthew G. Knepley PetscReal D = 1.0; 1672*29d8e7bcSMatthew G. Knepley PetscInt n, i; 1673*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1674*29d8e7bcSMatthew G. Knepley 1675*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1676*29d8e7bcSMatthew G. Knepley /* TODO: Assumes simplices */ 1677*29d8e7bcSMatthew G. Knepley ierr = DMPlexGetDimension(sp->dm, &n);CHKERRQ(ierr); 1678*29d8e7bcSMatthew G. Knepley for (i = 1; i <= n; ++i) { 1679*29d8e7bcSMatthew G. Knepley D *= ((PetscReal) (deg+i))/i; 1680*29d8e7bcSMatthew G. Knepley } 1681*29d8e7bcSMatthew G. Knepley *dim = (PetscInt) (D + 0.5); 1682*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1683*29d8e7bcSMatthew G. Knepley } 1684*29d8e7bcSMatthew G. Knepley 1685*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1686*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceInitialize_Lagrange" 1687*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 1688*29d8e7bcSMatthew G. Knepley { 1689*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1690*29d8e7bcSMatthew G. Knepley sp->ops->setfromoptions = NULL; 1691*29d8e7bcSMatthew G. Knepley sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 1692*29d8e7bcSMatthew G. Knepley sp->ops->view = NULL; 1693*29d8e7bcSMatthew G. Knepley sp->ops->destroy = NULL; 1694*29d8e7bcSMatthew G. Knepley sp->ops->getdimension = PetscDualSpaceGetDimension_Lagrange; 1695*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1696*29d8e7bcSMatthew G. Knepley } 1697*29d8e7bcSMatthew G. Knepley 1698*29d8e7bcSMatthew G. Knepley /*MC 1699*29d8e7bcSMatthew G. Knepley PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals 1700*29d8e7bcSMatthew G. Knepley 1701*29d8e7bcSMatthew G. Knepley Level: intermediate 1702*29d8e7bcSMatthew G. Knepley 1703*29d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 1704*29d8e7bcSMatthew G. Knepley M*/ 1705*29d8e7bcSMatthew G. Knepley 1706*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1707*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceCreate_Lagrange" 1708*29d8e7bcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 1709*29d8e7bcSMatthew G. Knepley { 1710*29d8e7bcSMatthew G. Knepley PetscDualSpace_Lag *lag; 1711*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1712*29d8e7bcSMatthew G. Knepley 1713*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1714*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1715*29d8e7bcSMatthew G. Knepley ierr = PetscNewLog(sp, PetscDualSpace_Lag, &lag);CHKERRQ(ierr); 1716*29d8e7bcSMatthew G. Knepley sp->data = lag; 1717*29d8e7bcSMatthew G. Knepley 1718*29d8e7bcSMatthew G. Knepley /* lag->n = 0; */ 1719*29d8e7bcSMatthew G. Knepley 1720*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr); 1721*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1722*29d8e7bcSMatthew G. Knepley } 1723*29d8e7bcSMatthew G. Knepley 1724*29d8e7bcSMatthew G. Knepley 1725*29d8e7bcSMatthew G. Knepley PetscInt PETSCFE_CLASSID = 0; 1726*29d8e7bcSMatthew G. Knepley 1727*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1728*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEView" 1729*29d8e7bcSMatthew G. Knepley /*@C 1730*29d8e7bcSMatthew G. Knepley PetscFEView - Views a PetscFE 1731*29d8e7bcSMatthew G. Knepley 1732*29d8e7bcSMatthew G. Knepley Collective on PetscFE 1733*29d8e7bcSMatthew G. Knepley 1734*29d8e7bcSMatthew G. Knepley Input Parameter: 1735*29d8e7bcSMatthew G. Knepley + sp - the PetscFE object to view 1736*29d8e7bcSMatthew G. Knepley - v - the viewer 1737*29d8e7bcSMatthew G. Knepley 1738*29d8e7bcSMatthew G. Knepley Level: developer 1739*29d8e7bcSMatthew G. Knepley 1740*29d8e7bcSMatthew G. Knepley .seealso PetscFEDestroy() 1741*29d8e7bcSMatthew G. Knepley @*/ 1742*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEView(PetscFE fem, PetscViewer v) 1743*29d8e7bcSMatthew G. Knepley { 1744*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1745*29d8e7bcSMatthew G. Knepley 1746*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1747*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1748*29d8e7bcSMatthew G. Knepley if (!v) { 1749*29d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &v);CHKERRQ(ierr); 1750*29d8e7bcSMatthew G. Knepley } 1751*29d8e7bcSMatthew G. Knepley if (fem->ops->view) { 1752*29d8e7bcSMatthew G. Knepley ierr = (*fem->ops->view)(fem, v);CHKERRQ(ierr); 1753*29d8e7bcSMatthew G. Knepley } 1754*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1755*29d8e7bcSMatthew G. Knepley } 1756*29d8e7bcSMatthew G. Knepley 1757*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1758*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEDestroy" 1759*29d8e7bcSMatthew G. Knepley /*@ 1760*29d8e7bcSMatthew G. Knepley PetscFEDestroy - Destroys a PetscFE object 1761*29d8e7bcSMatthew G. Knepley 1762*29d8e7bcSMatthew G. Knepley Collective on PetscFE 1763*29d8e7bcSMatthew G. Knepley 1764*29d8e7bcSMatthew G. Knepley Input Parameter: 1765*29d8e7bcSMatthew G. Knepley . fem - the PetscFE object to destroy 1766*29d8e7bcSMatthew G. Knepley 1767*29d8e7bcSMatthew G. Knepley Level: developer 1768*29d8e7bcSMatthew G. Knepley 1769*29d8e7bcSMatthew G. Knepley .seealso PetscFEView() 1770*29d8e7bcSMatthew G. Knepley @*/ 1771*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEDestroy(PetscFE *fem) 1772*29d8e7bcSMatthew G. Knepley { 1773*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1774*29d8e7bcSMatthew G. Knepley 1775*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1776*29d8e7bcSMatthew G. Knepley if (!*fem) PetscFunctionReturn(0); 1777*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1); 1778*29d8e7bcSMatthew G. Knepley 1779*29d8e7bcSMatthew G. Knepley if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);} 1780*29d8e7bcSMatthew G. Knepley ((PetscObject) (*fem))->refct = 0; 1781*29d8e7bcSMatthew G. Knepley /* if memory was published with AMS then destroy it */ 1782*29d8e7bcSMatthew G. Knepley ierr = PetscObjectAMSViewOff((PetscObject) *fem);CHKERRQ(ierr); 1783*29d8e7bcSMatthew G. Knepley 1784*29d8e7bcSMatthew G. Knepley ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr); 1785*29d8e7bcSMatthew G. Knepley ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr); 1786*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1787*29d8e7bcSMatthew G. Knepley } 1788*29d8e7bcSMatthew G. Knepley 1789*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1790*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFECreate" 1791*29d8e7bcSMatthew G. Knepley /*@ 1792*29d8e7bcSMatthew G. Knepley PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType(). 1793*29d8e7bcSMatthew G. Knepley 1794*29d8e7bcSMatthew G. Knepley Collective on MPI_Comm 1795*29d8e7bcSMatthew G. Knepley 1796*29d8e7bcSMatthew G. Knepley Input Parameter: 1797*29d8e7bcSMatthew G. Knepley . comm - The communicator for the PetscFE object 1798*29d8e7bcSMatthew G. Knepley 1799*29d8e7bcSMatthew G. Knepley Output Parameter: 1800*29d8e7bcSMatthew G. Knepley . fem - The PetscFE object 1801*29d8e7bcSMatthew G. Knepley 1802*29d8e7bcSMatthew G. Knepley Level: beginner 1803*29d8e7bcSMatthew G. Knepley 1804*29d8e7bcSMatthew G. Knepley .seealso: PetscFESetType(), PETSCFEGALERKIN 1805*29d8e7bcSMatthew G. Knepley @*/ 1806*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem) 1807*29d8e7bcSMatthew G. Knepley { 1808*29d8e7bcSMatthew G. Knepley PetscFE f; 1809*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1810*29d8e7bcSMatthew G. Knepley 1811*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1812*29d8e7bcSMatthew G. Knepley PetscValidPointer(fem, 2); 1813*29d8e7bcSMatthew G. Knepley *fem = NULL; 1814*29d8e7bcSMatthew G. Knepley #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1815*29d8e7bcSMatthew G. Knepley ierr = PetscFEInitializePackage();CHKERRQ(ierr); 1816*29d8e7bcSMatthew G. Knepley #endif 1817*29d8e7bcSMatthew G. Knepley 1818*29d8e7bcSMatthew G. Knepley ierr = PetscHeaderCreate(f, _p_PetscFE, struct _PetscFEOps, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr); 1819*29d8e7bcSMatthew G. Knepley ierr = PetscMemzero(f->ops, sizeof(struct _PetscFEOps));CHKERRQ(ierr); 1820*29d8e7bcSMatthew G. Knepley 1821*29d8e7bcSMatthew G. Knepley f->basisSpace = NULL; 1822*29d8e7bcSMatthew G. Knepley f->dualSpace = NULL; 1823*29d8e7bcSMatthew G. Knepley 1824*29d8e7bcSMatthew G. Knepley *fem = f; 1825*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1826*29d8e7bcSMatthew G. Knepley } 1827*29d8e7bcSMatthew G. Knepley 1828*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1829*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetDimension" 1830*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 1831*29d8e7bcSMatthew G. Knepley { 1832*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1833*29d8e7bcSMatthew G. Knepley 1834*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1835*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1836*29d8e7bcSMatthew G. Knepley PetscValidPointer(dim, 2); 1837*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceGetDimension(fem->basisSpace, dim);CHKERRQ(ierr); 1838*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1839*29d8e7bcSMatthew G. Knepley } 1840*29d8e7bcSMatthew G. Knepley 1841*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1842*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetNumComponents" 1843*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp) 1844*29d8e7bcSMatthew G. Knepley { 1845*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1846*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1847*29d8e7bcSMatthew G. Knepley PetscValidPointer(comp, 2); 1848*29d8e7bcSMatthew G. Knepley *comp = 1; 1849*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1850*29d8e7bcSMatthew G. Knepley } 1851*29d8e7bcSMatthew G. Knepley 1852*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1853*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetBasisSpace" 1854*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp) 1855*29d8e7bcSMatthew G. Knepley { 1856*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1857*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1858*29d8e7bcSMatthew G. Knepley PetscValidPointer(sp, 2); 1859*29d8e7bcSMatthew G. Knepley *sp = fem->basisSpace; 1860*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1861*29d8e7bcSMatthew G. Knepley } 1862*29d8e7bcSMatthew G. Knepley 1863*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1864*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFESetBasisSpace" 1865*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp) 1866*29d8e7bcSMatthew G. Knepley { 1867*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1868*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1869*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2); 1870*29d8e7bcSMatthew G. Knepley fem->basisSpace = sp; 1871*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1872*29d8e7bcSMatthew G. Knepley } 1873*29d8e7bcSMatthew G. Knepley 1874*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1875*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetDualSpace" 1876*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp) 1877*29d8e7bcSMatthew G. Knepley { 1878*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1879*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1880*29d8e7bcSMatthew G. Knepley PetscValidPointer(sp, 2); 1881*29d8e7bcSMatthew G. Knepley *sp = fem->dualSpace; 1882*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1883*29d8e7bcSMatthew G. Knepley } 1884*29d8e7bcSMatthew G. Knepley 1885*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1886*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFESetDualSpace" 1887*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp) 1888*29d8e7bcSMatthew G. Knepley { 1889*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1890*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1891*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 1892*29d8e7bcSMatthew G. Knepley fem->dualSpace = sp; 1893*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1894*29d8e7bcSMatthew G. Knepley } 1895*29d8e7bcSMatthew G. Knepley 1896*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1897*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetTabulation" 1898*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 1899*29d8e7bcSMatthew G. Knepley { 1900*29d8e7bcSMatthew G. Knepley DM dm; 1901*29d8e7bcSMatthew G. Knepley PetscReal *tmpB, *invV; 1902*29d8e7bcSMatthew G. Knepley PetscInt pdim; /* Dimension of FE space P */ 1903*29d8e7bcSMatthew G. Knepley PetscInt dim; /* Spatial dimension */ 1904*29d8e7bcSMatthew G. Knepley PetscInt p, j, k; 1905*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1906*29d8e7bcSMatthew G. Knepley 1907*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1908*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1909*29d8e7bcSMatthew G. Knepley PetscValidPointer(points, 3); 1910*29d8e7bcSMatthew G. Knepley if (B) PetscValidPointer(B, 4); 1911*29d8e7bcSMatthew G. Knepley if (D) PetscValidPointer(D, 5); 1912*29d8e7bcSMatthew G. Knepley if (H) PetscValidPointer(H, 6); 1913*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 1914*29d8e7bcSMatthew G. Knepley 1915*29d8e7bcSMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1916*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceGetDimension(fem->basisSpace, &pdim);CHKERRQ(ierr); 1917*29d8e7bcSMatthew G. Knepley /* 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*29d8e7bcSMatthew G. Knepley 1919*29d8e7bcSMatthew G. Knepley if (pdim != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Code only works for linear elements right now"); 1920*29d8e7bcSMatthew G. Knepley 1921*29d8e7bcSMatthew G. Knepley if (B) { 1922*29d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, npoints*pdim, PETSC_REAL, B);CHKERRQ(ierr); 1923*29d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr); 1924*29d8e7bcSMatthew G. Knepley } 1925*29d8e7bcSMatthew G. Knepley if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*dim, PETSC_REAL, D);CHKERRQ(ierr);} 1926*29d8e7bcSMatthew G. Knepley if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*dim*dim, PETSC_REAL, H);CHKERRQ(ierr);} 1927*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr); 1928*29d8e7bcSMatthew G. Knepley 1929*29d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr); 1930*29d8e7bcSMatthew G. Knepley for (j = 0; j < pdim; ++j) { 1931*29d8e7bcSMatthew G. Knepley PetscReal *Bf; 1932*29d8e7bcSMatthew G. Knepley PetscQuadrature f; 1933*29d8e7bcSMatthew G. Knepley PetscInt q; 1934*29d8e7bcSMatthew G. Knepley 1935*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f); 1936*29d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr); 1937*29d8e7bcSMatthew G. Knepley ierr = PetscSpaceEvaluate(fem->basisSpace, f.numQuadPoints, f.quadPoints, Bf, NULL, NULL);CHKERRQ(ierr); 1938*29d8e7bcSMatthew G. Knepley for (k = 0; k < pdim; ++k) { 1939*29d8e7bcSMatthew G. Knepley /* n_j \cdot \phi_k */ 1940*29d8e7bcSMatthew G. Knepley invV[j*pdim+k] = 0.0; 1941*29d8e7bcSMatthew G. Knepley for (q = 0; q < f.numQuadPoints; ++q) { 1942*29d8e7bcSMatthew G. Knepley invV[j*pdim+k] += Bf[q*pdim+k]*f.quadWeights[q]; 1943*29d8e7bcSMatthew G. Knepley } 1944*29d8e7bcSMatthew G. Knepley } 1945*29d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr); 1946*29d8e7bcSMatthew G. Knepley } 1947*29d8e7bcSMatthew G. Knepley { 1948*29d8e7bcSMatthew G. Knepley PetscReal *work; 1949*29d8e7bcSMatthew G. Knepley PetscBLASInt *pivots; 1950*29d8e7bcSMatthew G. Knepley PetscBLASInt n = pdim, info; 1951*29d8e7bcSMatthew G. Knepley 1952*29d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr); 1953*29d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr); 1954*29d8e7bcSMatthew G. Knepley PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invV, &n, pivots, &info)); 1955*29d8e7bcSMatthew G. Knepley PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invV, &n, pivots, work, &n, &info)); 1956*29d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr); 1957*29d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr); 1958*29d8e7bcSMatthew G. Knepley } 1959*29d8e7bcSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 1960*29d8e7bcSMatthew G. Knepley if (B) { 1961*29d8e7bcSMatthew G. Knepley /* Multiply by V^{-1} (pdim x pdim) */ 1962*29d8e7bcSMatthew G. Knepley for (j = 0; j < pdim; ++j) { 1963*29d8e7bcSMatthew G. Knepley const PetscInt i = p*pdim + j; 1964*29d8e7bcSMatthew G. Knepley 1965*29d8e7bcSMatthew G. Knepley (*B)[i] = 0.0; 1966*29d8e7bcSMatthew G. Knepley for (k = 0; k < pdim; ++k) { 1967*29d8e7bcSMatthew G. Knepley (*B)[i] += invV[k*pdim+j] * tmpB[p*pdim + k]; 1968*29d8e7bcSMatthew G. Knepley } 1969*29d8e7bcSMatthew G. Knepley } 1970*29d8e7bcSMatthew G. Knepley } 1971*29d8e7bcSMatthew G. Knepley } 1972*29d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr); 1973*29d8e7bcSMatthew G. Knepley if (B) { 1974*29d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr); 1975*29d8e7bcSMatthew G. Knepley } 1976*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1977*29d8e7bcSMatthew G. Knepley } 1978*29d8e7bcSMatthew G. Knepley 1979*29d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1980*29d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFERestoreTabulation" 1981*29d8e7bcSMatthew G. Knepley PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **D2) 1982*29d8e7bcSMatthew G. Knepley { 1983*29d8e7bcSMatthew G. Knepley DM dm; 1984*29d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 1985*29d8e7bcSMatthew G. Knepley 1986*29d8e7bcSMatthew G. Knepley PetscFunctionBegin; 1987*29d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1988*29d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 1989*29d8e7bcSMatthew G. Knepley if (B) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, B);CHKERRQ(ierr);} 1990*29d8e7bcSMatthew G. Knepley if (D) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D);CHKERRQ(ierr);} 1991*29d8e7bcSMatthew G. Knepley if (D2) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D2);CHKERRQ(ierr);} 1992*29d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 1993*29d8e7bcSMatthew G. Knepley } 1994