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 83903c00fSMatthew G. Knepley #include <petscdt.h> /*I "petscdt.h" I*/ /*I "petscfe.h" I*/ 937045ce4SJed Brown #include <petscblaslapack.h> 10194825f6SJed Brown #include <petsc-private/petscimpl.h> 11665c2dedSJed Brown #include <petscviewer.h> 1259804f93SMatthew G. Knepley #include <petscdmplex.h> 1359804f93SMatthew G. Knepley #include <petscdmshell.h> 1437045ce4SJed Brown 1537045ce4SJed Brown #undef __FUNCT__ 1637045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval" 1737045ce4SJed Brown /*@ 1837045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 1937045ce4SJed Brown 2037045ce4SJed Brown Not Collective 2137045ce4SJed Brown 2237045ce4SJed Brown Input Arguments: 2337045ce4SJed Brown + npoints - number of spatial points to evaluate at 2437045ce4SJed Brown . points - array of locations to evaluate at 2537045ce4SJed Brown . ndegree - number of basis degrees to evaluate 2637045ce4SJed Brown - degrees - sorted array of degrees to evaluate 2737045ce4SJed Brown 2837045ce4SJed Brown Output Arguments: 290298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 300298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 310298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 3237045ce4SJed Brown 3337045ce4SJed Brown Level: intermediate 3437045ce4SJed Brown 3537045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 3637045ce4SJed Brown @*/ 3737045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 3837045ce4SJed Brown { 3937045ce4SJed Brown PetscInt i,maxdegree; 4037045ce4SJed Brown 4137045ce4SJed Brown PetscFunctionBegin; 4237045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 4337045ce4SJed Brown maxdegree = degrees[ndegree-1]; 4437045ce4SJed Brown for (i=0; i<npoints; i++) { 4537045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 4637045ce4SJed Brown PetscInt j,k; 4737045ce4SJed Brown x = points[i]; 4837045ce4SJed Brown pm2 = 0; 4937045ce4SJed Brown pm1 = 1; 5037045ce4SJed Brown pd2 = 0; 5137045ce4SJed Brown pd1 = 0; 5237045ce4SJed Brown pdd2 = 0; 5337045ce4SJed Brown pdd1 = 0; 5437045ce4SJed Brown k = 0; 5537045ce4SJed Brown if (degrees[k] == 0) { 5637045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 5737045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 5837045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 5937045ce4SJed Brown k++; 6037045ce4SJed Brown } 6137045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 6237045ce4SJed Brown PetscReal p,d,dd; 6337045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 6437045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 6537045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 6637045ce4SJed Brown pm2 = pm1; 6737045ce4SJed Brown pm1 = p; 6837045ce4SJed Brown pd2 = pd1; 6937045ce4SJed Brown pd1 = d; 7037045ce4SJed Brown pdd2 = pdd1; 7137045ce4SJed Brown pdd1 = dd; 7237045ce4SJed Brown if (degrees[k] == j) { 7337045ce4SJed Brown if (B) B[i*ndegree+k] = p; 7437045ce4SJed Brown if (D) D[i*ndegree+k] = d; 7537045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 7637045ce4SJed Brown } 7737045ce4SJed Brown } 7837045ce4SJed Brown } 7937045ce4SJed Brown PetscFunctionReturn(0); 8037045ce4SJed Brown } 8137045ce4SJed Brown 8237045ce4SJed Brown #undef __FUNCT__ 8337045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature" 8437045ce4SJed Brown /*@ 8537045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 8637045ce4SJed Brown 8737045ce4SJed Brown Not Collective 8837045ce4SJed Brown 8937045ce4SJed Brown Input Arguments: 9037045ce4SJed Brown + npoints - number of points 9137045ce4SJed Brown . a - left end of interval (often-1) 9237045ce4SJed Brown - b - right end of interval (often +1) 9337045ce4SJed Brown 9437045ce4SJed Brown Output Arguments: 9537045ce4SJed Brown + x - quadrature points 9637045ce4SJed Brown - w - quadrature weights 9737045ce4SJed Brown 9837045ce4SJed Brown Level: intermediate 9937045ce4SJed Brown 10037045ce4SJed Brown References: 10137045ce4SJed Brown Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 10237045ce4SJed Brown 10337045ce4SJed Brown .seealso: PetscDTLegendreEval() 10437045ce4SJed Brown @*/ 10537045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 10637045ce4SJed Brown { 10737045ce4SJed Brown PetscErrorCode ierr; 10837045ce4SJed Brown PetscInt i; 10937045ce4SJed Brown PetscReal *work; 11037045ce4SJed Brown PetscScalar *Z; 11137045ce4SJed Brown PetscBLASInt N,LDZ,info; 11237045ce4SJed Brown 11337045ce4SJed Brown PetscFunctionBegin; 11437045ce4SJed Brown /* Set up the Golub-Welsch system */ 11537045ce4SJed Brown for (i=0; i<npoints; i++) { 11637045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 11737045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 11837045ce4SJed Brown } 11937045ce4SJed Brown ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 12037045ce4SJed Brown ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr); 121c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 12237045ce4SJed Brown LDZ = N; 12337045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1248b83055fSJed Brown PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 12537045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 1261c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 12737045ce4SJed Brown 12837045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 12937045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 13037045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 13137045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 1320d644c17SKarl Rupp 13337045ce4SJed Brown w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints])); 13437045ce4SJed Brown } 13537045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 13637045ce4SJed Brown PetscFunctionReturn(0); 13737045ce4SJed Brown } 138194825f6SJed Brown 139194825f6SJed Brown #undef __FUNCT__ 140494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal" 141494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 142494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 143494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial) 144494e7359SMatthew G. Knepley { 145494e7359SMatthew G. Knepley PetscReal f = 1.0; 146494e7359SMatthew G. Knepley PetscInt i; 147494e7359SMatthew G. Knepley 148494e7359SMatthew G. Knepley PetscFunctionBegin; 149494e7359SMatthew G. Knepley for (i = 1; i < n+1; ++i) f *= i; 150494e7359SMatthew G. Knepley *factorial = f; 151494e7359SMatthew G. Knepley PetscFunctionReturn(0); 152494e7359SMatthew G. Knepley } 153494e7359SMatthew G. Knepley 154494e7359SMatthew G. Knepley #undef __FUNCT__ 155494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi" 156494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 157494e7359SMatthew G. Knepley Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 158494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 159494e7359SMatthew G. Knepley { 160494e7359SMatthew G. Knepley PetscReal apb, pn1, pn2; 161494e7359SMatthew G. Knepley PetscInt k; 162494e7359SMatthew G. Knepley 163494e7359SMatthew G. Knepley PetscFunctionBegin; 164494e7359SMatthew G. Knepley if (!n) {*P = 1.0; PetscFunctionReturn(0);} 165494e7359SMatthew G. Knepley if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);} 166494e7359SMatthew G. Knepley apb = a + b; 167494e7359SMatthew G. Knepley pn2 = 1.0; 168494e7359SMatthew G. Knepley pn1 = 0.5 * (a - b + (apb + 2.0) * x); 169494e7359SMatthew G. Knepley *P = 0.0; 170494e7359SMatthew G. Knepley for (k = 2; k < n+1; ++k) { 171494e7359SMatthew G. Knepley PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0); 172494e7359SMatthew G. Knepley PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b); 173494e7359SMatthew G. Knepley PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb); 174494e7359SMatthew G. Knepley PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb); 175494e7359SMatthew G. Knepley 176494e7359SMatthew G. Knepley a2 = a2 / a1; 177494e7359SMatthew G. Knepley a3 = a3 / a1; 178494e7359SMatthew G. Knepley a4 = a4 / a1; 179494e7359SMatthew G. Knepley *P = (a2 + a3 * x) * pn1 - a4 * pn2; 180494e7359SMatthew G. Knepley pn2 = pn1; 181494e7359SMatthew G. Knepley pn1 = *P; 182494e7359SMatthew G. Knepley } 183494e7359SMatthew G. Knepley PetscFunctionReturn(0); 184494e7359SMatthew G. Knepley } 185494e7359SMatthew G. Knepley 186494e7359SMatthew G. Knepley #undef __FUNCT__ 187494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative" 188494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 189494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 190494e7359SMatthew G. Knepley { 191494e7359SMatthew G. Knepley PetscReal nP; 192494e7359SMatthew G. Knepley PetscErrorCode ierr; 193494e7359SMatthew G. Knepley 194494e7359SMatthew G. Knepley PetscFunctionBegin; 195494e7359SMatthew G. Knepley if (!n) {*P = 0.0; PetscFunctionReturn(0);} 196494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr); 197494e7359SMatthew G. Knepley *P = 0.5 * (a + b + n + 1) * nP; 198494e7359SMatthew G. Knepley PetscFunctionReturn(0); 199494e7359SMatthew G. Knepley } 200494e7359SMatthew G. Knepley 201494e7359SMatthew G. Knepley #undef __FUNCT__ 202494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal" 203494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 204494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 205494e7359SMatthew G. Knepley { 206494e7359SMatthew G. Knepley PetscFunctionBegin; 207494e7359SMatthew G. Knepley *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 208494e7359SMatthew G. Knepley *eta = y; 209494e7359SMatthew G. Knepley PetscFunctionReturn(0); 210494e7359SMatthew G. Knepley } 211494e7359SMatthew G. Knepley 212494e7359SMatthew G. Knepley #undef __FUNCT__ 213494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal" 214494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 215494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 216494e7359SMatthew G. Knepley { 217494e7359SMatthew G. Knepley PetscFunctionBegin; 218494e7359SMatthew G. Knepley *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 219494e7359SMatthew G. Knepley *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 220494e7359SMatthew G. Knepley *zeta = z; 221494e7359SMatthew G. Knepley PetscFunctionReturn(0); 222494e7359SMatthew G. Knepley } 223494e7359SMatthew G. Knepley 224494e7359SMatthew G. Knepley #undef __FUNCT__ 225494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal" 226494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 227494e7359SMatthew G. Knepley { 228494e7359SMatthew G. Knepley PetscInt maxIter = 100; 229494e7359SMatthew G. Knepley PetscReal eps = 1.0e-8; 230a8291ba1SSatish Balay PetscReal a1, a2, a3, a4, a5, a6; 231494e7359SMatthew G. Knepley PetscInt k; 232494e7359SMatthew G. Knepley PetscErrorCode ierr; 233494e7359SMatthew G. Knepley 234494e7359SMatthew G. Knepley PetscFunctionBegin; 235a8291ba1SSatish Balay 236a8291ba1SSatish Balay a1 = pow(2, a+b+1); 237a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA) 238a8291ba1SSatish Balay a2 = tgamma(a + npoints + 1); 239a8291ba1SSatish Balay a3 = tgamma(b + npoints + 1); 240a8291ba1SSatish Balay a4 = tgamma(a + b + npoints + 1); 241a8291ba1SSatish Balay #else 242a8291ba1SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 243a8291ba1SSatish Balay #endif 244a8291ba1SSatish Balay 245494e7359SMatthew G. Knepley ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr); 246494e7359SMatthew G. Knepley a6 = a1 * a2 * a3 / a4 / a5; 247494e7359SMatthew G. Knepley /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 248494e7359SMatthew G. Knepley Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 249494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 2507f1c68b3SMatthew G. Knepley PetscReal r = -cos((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP; 251494e7359SMatthew G. Knepley PetscInt j; 252494e7359SMatthew G. Knepley 253494e7359SMatthew G. Knepley if (k > 0) r = 0.5 * (r + x[k-1]); 254494e7359SMatthew G. Knepley for (j = 0; j < maxIter; ++j) { 255494e7359SMatthew G. Knepley PetscReal s = 0.0, delta, f, fp; 256494e7359SMatthew G. Knepley PetscInt i; 257494e7359SMatthew G. Knepley 258494e7359SMatthew G. Knepley for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 259494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 260494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); 261494e7359SMatthew G. Knepley delta = f / (fp - f * s); 262494e7359SMatthew G. Knepley r = r - delta; 263494e7359SMatthew G. Knepley if (fabs(delta) < eps) break; 264494e7359SMatthew G. Knepley } 265494e7359SMatthew G. Knepley x[k] = r; 266494e7359SMatthew G. Knepley ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); 267494e7359SMatthew G. Knepley w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 268494e7359SMatthew G. Knepley } 269494e7359SMatthew G. Knepley PetscFunctionReturn(0); 270494e7359SMatthew G. Knepley } 271494e7359SMatthew G. Knepley 272494e7359SMatthew G. Knepley #undef __FUNCT__ 273494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature" 274fd9d31fbSMatthew G. Knepley /*@C 275494e7359SMatthew G. Knepley PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex 276494e7359SMatthew G. Knepley 277494e7359SMatthew G. Knepley Not Collective 278494e7359SMatthew G. Knepley 279494e7359SMatthew G. Knepley Input Arguments: 280494e7359SMatthew G. Knepley + dim - The simplex dimension 281494e7359SMatthew G. Knepley . npoints - number of points 282494e7359SMatthew G. Knepley . a - left end of interval (often-1) 283494e7359SMatthew G. Knepley - b - right end of interval (often +1) 284494e7359SMatthew G. Knepley 285494e7359SMatthew G. Knepley Output Arguments: 286494e7359SMatthew G. Knepley + points - quadrature points 287494e7359SMatthew G. Knepley - weights - quadrature weights 288494e7359SMatthew G. Knepley 289494e7359SMatthew G. Knepley Level: intermediate 290494e7359SMatthew G. Knepley 291494e7359SMatthew G. Knepley References: 292494e7359SMatthew G. Knepley Karniadakis and Sherwin. 293494e7359SMatthew G. Knepley FIAT 294494e7359SMatthew G. Knepley 295494e7359SMatthew G. Knepley .seealso: PetscDTGaussQuadrature() 296494e7359SMatthew G. Knepley @*/ 297494e7359SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscReal *points[], PetscReal *weights[]) 298494e7359SMatthew G. Knepley { 299494e7359SMatthew G. Knepley PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 300494e7359SMatthew G. Knepley PetscInt i, j, k; 301494e7359SMatthew G. Knepley PetscErrorCode ierr; 302494e7359SMatthew G. Knepley 303494e7359SMatthew G. Knepley PetscFunctionBegin; 304494e7359SMatthew G. Knepley if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 305494e7359SMatthew G. Knepley switch (dim) { 306494e7359SMatthew G. Knepley case 1: 307494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints * sizeof(PetscReal), &x);CHKERRQ(ierr); 308494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 309494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, w);CHKERRQ(ierr); 310494e7359SMatthew G. Knepley break; 311494e7359SMatthew G. Knepley case 2: 312494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints*npoints*2 * sizeof(PetscReal), &x);CHKERRQ(ierr); 313494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints*npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 314494e7359SMatthew G. Knepley ierr = PetscMalloc4(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy);CHKERRQ(ierr); 315494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 316494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 317494e7359SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 318494e7359SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 319494e7359SMatthew G. Knepley ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 320494e7359SMatthew G. Knepley w[i*npoints+j] = 0.5 * wx[i] * wy[j]; 321494e7359SMatthew G. Knepley } 322494e7359SMatthew G. Knepley } 323494e7359SMatthew G. Knepley ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 324494e7359SMatthew G. Knepley break; 325494e7359SMatthew G. Knepley case 3: 326494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints*npoints*3 * sizeof(PetscReal), &x);CHKERRQ(ierr); 327494e7359SMatthew G. Knepley ierr = PetscMalloc(npoints*npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); 328494e7359SMatthew G. Knepley ierr = PetscMalloc6(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy,npoints,PetscReal,&pz,npoints,PetscReal,&wz);CHKERRQ(ierr); 329494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); 330494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); 331494e7359SMatthew G. Knepley ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); 332494e7359SMatthew G. Knepley for (i = 0; i < npoints; ++i) { 333494e7359SMatthew G. Knepley for (j = 0; j < npoints; ++j) { 334494e7359SMatthew G. Knepley for (k = 0; k < npoints; ++k) { 335494e7359SMatthew 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); 336494e7359SMatthew G. Knepley w[(i*npoints+j)*npoints+k] = 0.125 * wx[i] * wy[j] * wz[k]; 337494e7359SMatthew G. Knepley } 338494e7359SMatthew G. Knepley } 339494e7359SMatthew G. Knepley } 340494e7359SMatthew G. Knepley ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 341494e7359SMatthew G. Knepley break; 342494e7359SMatthew G. Knepley default: 343494e7359SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 344494e7359SMatthew G. Knepley } 345494e7359SMatthew G. Knepley if (points) *points = x; 346494e7359SMatthew G. Knepley if (weights) *weights = w; 347494e7359SMatthew G. Knepley PetscFunctionReturn(0); 348494e7359SMatthew G. Knepley } 349494e7359SMatthew G. Knepley 350494e7359SMatthew G. Knepley #undef __FUNCT__ 351194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR" 352194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 353194825f6SJed Brown * A in column-major format 354194825f6SJed Brown * Ainv in row-major format 355194825f6SJed Brown * tau has length m 356194825f6SJed Brown * worksize must be >= max(1,n) 357194825f6SJed Brown */ 358194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 359194825f6SJed Brown { 360194825f6SJed Brown PetscErrorCode ierr; 361194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 362194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 363194825f6SJed Brown 364194825f6SJed Brown PetscFunctionBegin; 365194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 366194825f6SJed Brown { 367194825f6SJed Brown PetscInt i,j; 368194825f6SJed Brown ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr); 369194825f6SJed Brown for (j=0; j<n; j++) { 370194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 371194825f6SJed Brown } 372194825f6SJed Brown mstride = m; 373194825f6SJed Brown } 374194825f6SJed Brown #else 375194825f6SJed Brown A = A_in; 376194825f6SJed Brown Ainv = Ainv_out; 377194825f6SJed Brown #endif 378194825f6SJed Brown 379194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 380194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 381194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 382194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 383194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 384194825f6SJed Brown LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info); 385194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 386194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 387194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 388194825f6SJed Brown 389194825f6SJed Brown /* Extract an explicit representation of Q */ 390194825f6SJed Brown Q = Ainv; 391194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 392194825f6SJed Brown K = N; /* full rank */ 393194825f6SJed Brown LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info); 394194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 395194825f6SJed Brown 396194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 397194825f6SJed Brown Alpha = 1.0; 398194825f6SJed Brown ldb = lda; 399194825f6SJed Brown BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb); 400194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 401194825f6SJed Brown 402194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 403194825f6SJed Brown { 404194825f6SJed Brown PetscInt i; 405194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 406194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 407194825f6SJed Brown } 408194825f6SJed Brown #endif 409194825f6SJed Brown PetscFunctionReturn(0); 410194825f6SJed Brown } 411194825f6SJed Brown 412194825f6SJed Brown #undef __FUNCT__ 413194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate" 414194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 415194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 416194825f6SJed Brown { 417194825f6SJed Brown PetscErrorCode ierr; 418194825f6SJed Brown PetscReal *Bv; 419194825f6SJed Brown PetscInt i,j; 420194825f6SJed Brown 421194825f6SJed Brown PetscFunctionBegin; 422194825f6SJed Brown ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr); 423194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 424194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 425194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 426194825f6SJed Brown for (i=0; i<ninterval; i++) { 427194825f6SJed Brown for (j=0; j<ndegree; j++) { 428194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 429194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 430194825f6SJed Brown } 431194825f6SJed Brown } 432194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 433194825f6SJed Brown PetscFunctionReturn(0); 434194825f6SJed Brown } 435194825f6SJed Brown 436194825f6SJed Brown #undef __FUNCT__ 437194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly" 438194825f6SJed Brown /*@ 439194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 440194825f6SJed Brown 441194825f6SJed Brown Not Collective 442194825f6SJed Brown 443194825f6SJed Brown Input Arguments: 444194825f6SJed Brown + degree - degree of reconstruction polynomial 445194825f6SJed Brown . nsource - number of source intervals 446194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 447194825f6SJed Brown . ntarget - number of target intervals 448194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 449194825f6SJed Brown 450194825f6SJed Brown Output Arguments: 451194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 452194825f6SJed Brown 453194825f6SJed Brown Level: advanced 454194825f6SJed Brown 455194825f6SJed Brown .seealso: PetscDTLegendreEval() 456194825f6SJed Brown @*/ 457194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 458194825f6SJed Brown { 459194825f6SJed Brown PetscErrorCode ierr; 460194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 461194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 462194825f6SJed Brown PetscScalar *tau,*work; 463194825f6SJed Brown 464194825f6SJed Brown PetscFunctionBegin; 465194825f6SJed Brown PetscValidRealPointer(sourcex,3); 466194825f6SJed Brown PetscValidRealPointer(targetx,5); 467194825f6SJed Brown PetscValidRealPointer(R,6); 468194825f6SJed 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); 469194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 470194825f6SJed Brown for (i=0; i<nsource; i++) { 471194825f6SJed 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]); 472194825f6SJed Brown } 473194825f6SJed Brown for (i=0; i<ntarget; i++) { 474194825f6SJed 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]); 475194825f6SJed Brown } 476194825f6SJed Brown #endif 477194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 478194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 479194825f6SJed Brown center = (xmin + xmax)/2; 480194825f6SJed Brown hscale = (xmax - xmin)/2; 481194825f6SJed Brown worksize = nsource; 482194825f6SJed Brown ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr); 48382772646SJed Brown ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr); 484194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 485194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 486194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 487194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 488194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 489194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 490194825f6SJed Brown for (i=0; i<ntarget; i++) { 491194825f6SJed Brown PetscReal rowsum = 0; 492194825f6SJed Brown for (j=0; j<nsource; j++) { 493194825f6SJed Brown PetscReal sum = 0; 494194825f6SJed Brown for (k=0; k<degree+1; k++) { 495194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 496194825f6SJed Brown } 497194825f6SJed Brown R[i*nsource+j] = sum; 498194825f6SJed Brown rowsum += sum; 499194825f6SJed Brown } 500194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 501194825f6SJed Brown } 502194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 503194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 504194825f6SJed Brown PetscFunctionReturn(0); 505194825f6SJed Brown } 50629d8e7bcSMatthew G. Knepley 50729d8e7bcSMatthew G. Knepley /* Basis Jet Tabulation 50829d8e7bcSMatthew G. Knepley 50929d8e7bcSMatthew G. Knepley We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We 51029d8e7bcSMatthew G. Knepley follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can 51129d8e7bcSMatthew 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 51229d8e7bcSMatthew G. Knepley as a prime basis. 51329d8e7bcSMatthew G. Knepley 51429d8e7bcSMatthew G. Knepley \psi_i = \sum_k \alpha_{ki} \phi_k 51529d8e7bcSMatthew G. Knepley 51629d8e7bcSMatthew G. Knepley Our nodal basis is defined in terms of the dual basis $n_j$ 51729d8e7bcSMatthew G. Knepley 51829d8e7bcSMatthew G. Knepley n_j \cdot \psi_i = \delta_{ji} 51929d8e7bcSMatthew G. Knepley 52029d8e7bcSMatthew G. Knepley and we may act on the first equation to obtain 52129d8e7bcSMatthew G. Knepley 52229d8e7bcSMatthew G. Knepley n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k 52329d8e7bcSMatthew G. Knepley \delta_{ji} = \sum_k \alpha_{ki} V_{jk} 52429d8e7bcSMatthew G. Knepley I = V \alpha 52529d8e7bcSMatthew G. Knepley 52629d8e7bcSMatthew G. Knepley so the coefficients of the nodal basis in the prime basis are 52729d8e7bcSMatthew G. Knepley 52829d8e7bcSMatthew G. Knepley \alpha = V^{-1} 52929d8e7bcSMatthew G. Knepley 53029d8e7bcSMatthew G. Knepley We will define the dual basis vectors $n_j$ using a quadrature rule. 53129d8e7bcSMatthew G. Knepley 53229d8e7bcSMatthew G. Knepley Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials 53329d8e7bcSMatthew G. Knepley (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can 53429d8e7bcSMatthew G. Knepley be implemented exactly as in FIAT using functionals $L_j$. 53529d8e7bcSMatthew G. Knepley 53629d8e7bcSMatthew G. Knepley I will have to count the degrees correctly for the Legendre product when we are on simplices. 53729d8e7bcSMatthew G. Knepley 53829d8e7bcSMatthew G. Knepley We will have three objects: 53929d8e7bcSMatthew G. Knepley - Space, P: this just need point evaluation I think 54029d8e7bcSMatthew 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 54129d8e7bcSMatthew G. Knepley - FEM: This keeps {P, P', Q} 54229d8e7bcSMatthew G. Knepley */ 54329d8e7bcSMatthew G. Knepley #include <petsc-private/petscfeimpl.h> 54429d8e7bcSMatthew G. Knepley 54529d8e7bcSMatthew G. Knepley PetscInt PETSCSPACE_CLASSID = 0; 54629d8e7bcSMatthew G. Knepley 54729d8e7bcSMatthew G. Knepley PetscFunctionList PetscSpaceList = NULL; 54829d8e7bcSMatthew G. Knepley PetscBool PetscSpaceRegisterAllCalled = PETSC_FALSE; 54929d8e7bcSMatthew G. Knepley 55029d8e7bcSMatthew G. Knepley #undef __FUNCT__ 55129d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceRegister" 55229d8e7bcSMatthew G. Knepley /*@C 55329d8e7bcSMatthew G. Knepley PetscSpaceRegister - Adds a new PetscSpace implementation 55429d8e7bcSMatthew G. Knepley 55529d8e7bcSMatthew G. Knepley Not Collective 55629d8e7bcSMatthew G. Knepley 55729d8e7bcSMatthew G. Knepley Input Parameters: 55829d8e7bcSMatthew G. Knepley + name - The name of a new user-defined creation routine 55929d8e7bcSMatthew G. Knepley - create_func - The creation routine itself 56029d8e7bcSMatthew G. Knepley 56129d8e7bcSMatthew G. Knepley Notes: 56229d8e7bcSMatthew G. Knepley PetscSpaceRegister() may be called multiple times to add several user-defined PetscSpaces 56329d8e7bcSMatthew G. Knepley 56429d8e7bcSMatthew G. Knepley Sample usage: 56529d8e7bcSMatthew G. Knepley .vb 56629d8e7bcSMatthew G. Knepley PetscSpaceRegister("my_space", MyPetscSpaceCreate); 56729d8e7bcSMatthew G. Knepley .ve 56829d8e7bcSMatthew G. Knepley 56929d8e7bcSMatthew G. Knepley Then, your PetscSpace type can be chosen with the procedural interface via 57029d8e7bcSMatthew G. Knepley .vb 57129d8e7bcSMatthew G. Knepley PetscSpaceCreate(MPI_Comm, PetscSpace *); 57229d8e7bcSMatthew G. Knepley PetscSpaceSetType(PetscSpace, "my_space"); 57329d8e7bcSMatthew G. Knepley .ve 57429d8e7bcSMatthew G. Knepley or at runtime via the option 57529d8e7bcSMatthew G. Knepley .vb 57629d8e7bcSMatthew G. Knepley -petscspace_type my_space 57729d8e7bcSMatthew G. Knepley .ve 57829d8e7bcSMatthew G. Knepley 57929d8e7bcSMatthew G. Knepley Level: advanced 58029d8e7bcSMatthew G. Knepley 58129d8e7bcSMatthew G. Knepley .keywords: PetscSpace, register 58229d8e7bcSMatthew G. Knepley .seealso: PetscSpaceRegisterAll(), PetscSpaceRegisterDestroy() 58329d8e7bcSMatthew G. Knepley 58429d8e7bcSMatthew G. Knepley @*/ 58529d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscSpace)) 58629d8e7bcSMatthew G. Knepley { 58729d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 58829d8e7bcSMatthew G. Knepley 58929d8e7bcSMatthew G. Knepley PetscFunctionBegin; 59029d8e7bcSMatthew G. Knepley ierr = PetscFunctionListAdd(&PetscSpaceList, sname, function);CHKERRQ(ierr); 59129d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 59229d8e7bcSMatthew G. Knepley } 59329d8e7bcSMatthew G. Knepley 59429d8e7bcSMatthew G. Knepley #undef __FUNCT__ 59529d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetType" 59629d8e7bcSMatthew G. Knepley /*@C 59729d8e7bcSMatthew G. Knepley PetscSpaceSetType - Builds a particular PetscSpace 59829d8e7bcSMatthew G. Knepley 59929d8e7bcSMatthew G. Knepley Collective on PetscSpace 60029d8e7bcSMatthew G. Knepley 60129d8e7bcSMatthew G. Knepley Input Parameters: 60229d8e7bcSMatthew G. Knepley + sp - The PetscSpace object 60329d8e7bcSMatthew G. Knepley - name - The kind of space 60429d8e7bcSMatthew G. Knepley 60529d8e7bcSMatthew G. Knepley Options Database Key: 60629d8e7bcSMatthew G. Knepley . -petscspace_type <type> - Sets the PetscSpace type; use -help for a list of available types 60729d8e7bcSMatthew G. Knepley 60829d8e7bcSMatthew G. Knepley Level: intermediate 60929d8e7bcSMatthew G. Knepley 61029d8e7bcSMatthew G. Knepley .keywords: PetscSpace, set, type 61129d8e7bcSMatthew G. Knepley .seealso: PetscSpaceGetType(), PetscSpaceCreate() 61229d8e7bcSMatthew G. Knepley @*/ 61329d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetType(PetscSpace sp, PetscSpaceType name) 61429d8e7bcSMatthew G. Knepley { 61529d8e7bcSMatthew G. Knepley PetscErrorCode (*r)(PetscSpace); 61629d8e7bcSMatthew G. Knepley PetscBool match; 61729d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 61829d8e7bcSMatthew G. Knepley 61929d8e7bcSMatthew G. Knepley PetscFunctionBegin; 62029d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 62129d8e7bcSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 62229d8e7bcSMatthew G. Knepley if (match) PetscFunctionReturn(0); 62329d8e7bcSMatthew G. Knepley 62429d8e7bcSMatthew G. Knepley if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 62529d8e7bcSMatthew G. Knepley ierr = PetscFunctionListFind(PetscSpaceList, name, &r);CHKERRQ(ierr); 62629d8e7bcSMatthew G. Knepley if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSpace type: %s", name); 62729d8e7bcSMatthew G. Knepley 62829d8e7bcSMatthew G. Knepley if (sp->ops->destroy) { 62929d8e7bcSMatthew G. Knepley ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 63029d8e7bcSMatthew G. Knepley sp->ops->destroy = NULL; 63129d8e7bcSMatthew G. Knepley } 63229d8e7bcSMatthew G. Knepley ierr = (*r)(sp);CHKERRQ(ierr); 63329d8e7bcSMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 63429d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 63529d8e7bcSMatthew G. Knepley } 63629d8e7bcSMatthew G. Knepley 63729d8e7bcSMatthew G. Knepley #undef __FUNCT__ 63829d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceGetType" 63929d8e7bcSMatthew G. Knepley /*@C 64029d8e7bcSMatthew G. Knepley PetscSpaceGetType - Gets the PetscSpace type name (as a string) from the object. 64129d8e7bcSMatthew G. Knepley 64229d8e7bcSMatthew G. Knepley Not Collective 64329d8e7bcSMatthew G. Knepley 64429d8e7bcSMatthew G. Knepley Input Parameter: 64529d8e7bcSMatthew G. Knepley . dm - The PetscSpace 64629d8e7bcSMatthew G. Knepley 64729d8e7bcSMatthew G. Knepley Output Parameter: 64829d8e7bcSMatthew G. Knepley . name - The PetscSpace type name 64929d8e7bcSMatthew G. Knepley 65029d8e7bcSMatthew G. Knepley Level: intermediate 65129d8e7bcSMatthew G. Knepley 65229d8e7bcSMatthew G. Knepley .keywords: PetscSpace, get, type, name 65329d8e7bcSMatthew G. Knepley .seealso: PetscSpaceSetType(), PetscSpaceCreate() 65429d8e7bcSMatthew G. Knepley @*/ 65529d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceGetType(PetscSpace sp, PetscSpaceType *name) 65629d8e7bcSMatthew G. Knepley { 65729d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 65829d8e7bcSMatthew G. Knepley 65929d8e7bcSMatthew G. Knepley PetscFunctionBegin; 66029d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 66129d8e7bcSMatthew G. Knepley PetscValidCharPointer(name, 2); 66229d8e7bcSMatthew G. Knepley if (!PetscSpaceRegisterAllCalled) { 66329d8e7bcSMatthew G. Knepley ierr = PetscSpaceRegisterAll();CHKERRQ(ierr); 66429d8e7bcSMatthew G. Knepley } 66529d8e7bcSMatthew G. Knepley *name = ((PetscObject) sp)->type_name; 66629d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 66729d8e7bcSMatthew G. Knepley } 66829d8e7bcSMatthew G. Knepley 66929d8e7bcSMatthew G. Knepley #undef __FUNCT__ 67029d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceView" 67129d8e7bcSMatthew G. Knepley /*@C 67229d8e7bcSMatthew G. Knepley PetscSpaceView - Views a PetscSpace 67329d8e7bcSMatthew G. Knepley 67429d8e7bcSMatthew G. Knepley Collective on PetscSpace 67529d8e7bcSMatthew G. Knepley 67629d8e7bcSMatthew G. Knepley Input Parameter: 67729d8e7bcSMatthew G. Knepley + sp - the PetscSpace object to view 67829d8e7bcSMatthew G. Knepley - v - the viewer 67929d8e7bcSMatthew G. Knepley 68029d8e7bcSMatthew G. Knepley Level: developer 68129d8e7bcSMatthew G. Knepley 68229d8e7bcSMatthew G. Knepley .seealso PetscSpaceDestroy() 68329d8e7bcSMatthew G. Knepley @*/ 68429d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceView(PetscSpace sp, PetscViewer v) 68529d8e7bcSMatthew G. Knepley { 68629d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 68729d8e7bcSMatthew G. Knepley 68829d8e7bcSMatthew G. Knepley PetscFunctionBegin; 68929d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 69029d8e7bcSMatthew G. Knepley if (!v) { 69129d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr); 69229d8e7bcSMatthew G. Knepley } 69329d8e7bcSMatthew G. Knepley if (sp->ops->view) { 69429d8e7bcSMatthew G. Knepley ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr); 69529d8e7bcSMatthew G. Knepley } 69629d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 69729d8e7bcSMatthew G. Knepley } 69829d8e7bcSMatthew G. Knepley 69929d8e7bcSMatthew G. Knepley #undef __FUNCT__ 70029d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceViewFromOptions" 70129d8e7bcSMatthew G. Knepley /* 70229d8e7bcSMatthew G. Knepley PetscSpaceViewFromOptions - Processes command line options to determine if/how a PetscSpace is to be viewed. 70329d8e7bcSMatthew G. Knepley 70429d8e7bcSMatthew G. Knepley Collective on PetscSpace 70529d8e7bcSMatthew G. Knepley 70629d8e7bcSMatthew G. Knepley Input Parameters: 70729d8e7bcSMatthew G. Knepley + sp - the PetscSpace 70829d8e7bcSMatthew G. Knepley . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd' 70929d8e7bcSMatthew G. Knepley - optionname - option to activate viewing 71029d8e7bcSMatthew G. Knepley 71129d8e7bcSMatthew G. Knepley Level: intermediate 71229d8e7bcSMatthew G. Knepley 71329d8e7bcSMatthew G. Knepley .keywords: PetscSpace, view, options, database 71429d8e7bcSMatthew G. Knepley .seealso: VecViewFromOptions(), MatViewFromOptions() 71529d8e7bcSMatthew G. Knepley */ 71629d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceViewFromOptions(PetscSpace sp, const char prefix[], const char optionname[]) 71729d8e7bcSMatthew G. Knepley { 71829d8e7bcSMatthew G. Knepley PetscViewer viewer; 71929d8e7bcSMatthew G. Knepley PetscViewerFormat format; 72029d8e7bcSMatthew G. Knepley PetscBool flg; 72129d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 72229d8e7bcSMatthew G. Knepley 72329d8e7bcSMatthew G. Knepley PetscFunctionBegin; 72429d8e7bcSMatthew G. Knepley if (prefix) { 72529d8e7bcSMatthew G. Knepley ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 72629d8e7bcSMatthew G. Knepley } else { 72729d8e7bcSMatthew G. Knepley ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 72829d8e7bcSMatthew G. Knepley } 72929d8e7bcSMatthew G. Knepley if (flg) { 73029d8e7bcSMatthew G. Knepley ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 73129d8e7bcSMatthew G. Knepley ierr = PetscSpaceView(sp, viewer);CHKERRQ(ierr); 73229d8e7bcSMatthew G. Knepley ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 73329d8e7bcSMatthew G. Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 73429d8e7bcSMatthew G. Knepley } 73529d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 73629d8e7bcSMatthew G. Knepley } 73729d8e7bcSMatthew G. Knepley 73829d8e7bcSMatthew G. Knepley #undef __FUNCT__ 73929d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetFromOptions" 74029d8e7bcSMatthew G. Knepley /*@ 74129d8e7bcSMatthew G. Knepley PetscSpaceSetFromOptions - sets parameters in a PetscSpace from the options database 74229d8e7bcSMatthew G. Knepley 74329d8e7bcSMatthew G. Knepley Collective on PetscSpace 74429d8e7bcSMatthew G. Knepley 74529d8e7bcSMatthew G. Knepley Input Parameter: 74629d8e7bcSMatthew G. Knepley . sp - the PetscSpace object to set options for 74729d8e7bcSMatthew G. Knepley 74829d8e7bcSMatthew G. Knepley Options Database: 74929d8e7bcSMatthew G. Knepley . -petscspace_order the approximation order of the space 75029d8e7bcSMatthew G. Knepley 75129d8e7bcSMatthew G. Knepley Level: developer 75229d8e7bcSMatthew G. Knepley 75329d8e7bcSMatthew G. Knepley .seealso PetscSpaceView() 75429d8e7bcSMatthew G. Knepley @*/ 75529d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetFromOptions(PetscSpace sp) 75629d8e7bcSMatthew G. Knepley { 75729d8e7bcSMatthew G. Knepley const char *defaultType; 75829d8e7bcSMatthew G. Knepley char typename[256]; 75929d8e7bcSMatthew G. Knepley PetscBool flg; 76029d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 76129d8e7bcSMatthew G. Knepley 76229d8e7bcSMatthew G. Knepley PetscFunctionBegin; 76329d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 76429d8e7bcSMatthew G. Knepley if (!((PetscObject) sp)->type_name) { 76529d8e7bcSMatthew G. Knepley defaultType = PETSCSPACEPOLYNOMIAL; 76629d8e7bcSMatthew G. Knepley } else { 76729d8e7bcSMatthew G. Knepley defaultType = ((PetscObject) sp)->type_name; 76829d8e7bcSMatthew G. Knepley } 76929d8e7bcSMatthew G. Knepley if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 77029d8e7bcSMatthew G. Knepley 77129d8e7bcSMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 77229d8e7bcSMatthew G. Knepley ierr = PetscOptionsList("-petscspace_type", "Linear space", "PetscSpaceSetType", PetscSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr); 77329d8e7bcSMatthew G. Knepley if (flg) { 77429d8e7bcSMatthew G. Knepley ierr = PetscSpaceSetType(sp, typename);CHKERRQ(ierr); 77529d8e7bcSMatthew G. Knepley } else if (!((PetscObject) sp)->type_name) { 77629d8e7bcSMatthew G. Knepley ierr = PetscSpaceSetType(sp, defaultType);CHKERRQ(ierr); 77729d8e7bcSMatthew G. Knepley } 77829d8e7bcSMatthew G. Knepley ierr = PetscOptionsInt("-petscspace_order", "The approximation order", "PetscSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr); 77929d8e7bcSMatthew G. Knepley if (sp->ops->setfromoptions) { 78029d8e7bcSMatthew G. Knepley ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr); 78129d8e7bcSMatthew G. Knepley } 78229d8e7bcSMatthew G. Knepley /* process any options handlers added with PetscObjectAddOptionsHandler() */ 78329d8e7bcSMatthew G. Knepley ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr); 78429d8e7bcSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 78529d8e7bcSMatthew G. Knepley ierr = PetscSpaceViewFromOptions(sp, NULL, "-petscspace_view");CHKERRQ(ierr); 78629d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 78729d8e7bcSMatthew G. Knepley } 78829d8e7bcSMatthew G. Knepley 78929d8e7bcSMatthew G. Knepley #undef __FUNCT__ 79029d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetUp" 79129d8e7bcSMatthew G. Knepley /*@C 79229d8e7bcSMatthew G. Knepley PetscSpaceSetUp - Construct data structures for the PetscSpace 79329d8e7bcSMatthew G. Knepley 79429d8e7bcSMatthew G. Knepley Collective on PetscSpace 79529d8e7bcSMatthew G. Knepley 79629d8e7bcSMatthew G. Knepley Input Parameter: 79729d8e7bcSMatthew G. Knepley . sp - the PetscSpace object to setup 79829d8e7bcSMatthew G. Knepley 79929d8e7bcSMatthew G. Knepley Level: developer 80029d8e7bcSMatthew G. Knepley 80129d8e7bcSMatthew G. Knepley .seealso PetscSpaceView(), PetscSpaceDestroy() 80229d8e7bcSMatthew G. Knepley @*/ 80329d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetUp(PetscSpace sp) 80429d8e7bcSMatthew G. Knepley { 80529d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 80629d8e7bcSMatthew G. Knepley 80729d8e7bcSMatthew G. Knepley PetscFunctionBegin; 80829d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 80929d8e7bcSMatthew G. Knepley if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 81029d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 81129d8e7bcSMatthew G. Knepley } 81229d8e7bcSMatthew G. Knepley 81329d8e7bcSMatthew G. Knepley #undef __FUNCT__ 81429d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceDestroy" 81529d8e7bcSMatthew G. Knepley /*@ 81629d8e7bcSMatthew G. Knepley PetscSpaceDestroy - Destroys a PetscSpace object 81729d8e7bcSMatthew G. Knepley 81829d8e7bcSMatthew G. Knepley Collective on PetscSpace 81929d8e7bcSMatthew G. Knepley 82029d8e7bcSMatthew G. Knepley Input Parameter: 82129d8e7bcSMatthew G. Knepley . sp - the PetscSpace object to destroy 82229d8e7bcSMatthew G. Knepley 82329d8e7bcSMatthew G. Knepley Level: developer 82429d8e7bcSMatthew G. Knepley 82529d8e7bcSMatthew G. Knepley .seealso PetscSpaceView() 82629d8e7bcSMatthew G. Knepley @*/ 82729d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceDestroy(PetscSpace *sp) 82829d8e7bcSMatthew G. Knepley { 82929d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 83029d8e7bcSMatthew G. Knepley 83129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 83229d8e7bcSMatthew G. Knepley if (!*sp) PetscFunctionReturn(0); 83329d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific((*sp), PETSCSPACE_CLASSID, 1); 83429d8e7bcSMatthew G. Knepley 83529d8e7bcSMatthew G. Knepley if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 83629d8e7bcSMatthew G. Knepley ((PetscObject) (*sp))->refct = 0; 83729d8e7bcSMatthew G. Knepley /* if memory was published with AMS then destroy it */ 83829d8e7bcSMatthew G. Knepley ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr); 83929d8e7bcSMatthew G. Knepley 84029d8e7bcSMatthew G. Knepley ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr); 84129d8e7bcSMatthew G. Knepley 84229d8e7bcSMatthew G. Knepley ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr); 84329d8e7bcSMatthew G. Knepley ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 84429d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 84529d8e7bcSMatthew G. Knepley } 84629d8e7bcSMatthew G. Knepley 84729d8e7bcSMatthew G. Knepley #undef __FUNCT__ 84829d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceCreate" 84929d8e7bcSMatthew G. Knepley /*@ 85029d8e7bcSMatthew G. Knepley PetscSpaceCreate - Creates an empty PetscSpace object. The type can then be set with PetscSpaceSetType(). 85129d8e7bcSMatthew G. Knepley 85229d8e7bcSMatthew G. Knepley Collective on MPI_Comm 85329d8e7bcSMatthew G. Knepley 85429d8e7bcSMatthew G. Knepley Input Parameter: 85529d8e7bcSMatthew G. Knepley . comm - The communicator for the PetscSpace object 85629d8e7bcSMatthew G. Knepley 85729d8e7bcSMatthew G. Knepley Output Parameter: 85829d8e7bcSMatthew G. Knepley . sp - The PetscSpace object 85929d8e7bcSMatthew G. Knepley 86029d8e7bcSMatthew G. Knepley Level: beginner 86129d8e7bcSMatthew G. Knepley 86229d8e7bcSMatthew G. Knepley .seealso: PetscSpaceSetType(), PETSCSPACEPOLYNOMIAL 86329d8e7bcSMatthew G. Knepley @*/ 86429d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceCreate(MPI_Comm comm, PetscSpace *sp) 86529d8e7bcSMatthew G. Knepley { 86629d8e7bcSMatthew G. Knepley PetscSpace s; 86729d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 86829d8e7bcSMatthew G. Knepley 86929d8e7bcSMatthew G. Knepley PetscFunctionBegin; 87029d8e7bcSMatthew G. Knepley PetscValidPointer(sp, 2); 87129d8e7bcSMatthew G. Knepley *sp = NULL; 87229d8e7bcSMatthew G. Knepley #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 87329d8e7bcSMatthew G. Knepley ierr = PetscFEInitializePackage();CHKERRQ(ierr); 87429d8e7bcSMatthew G. Knepley #endif 87529d8e7bcSMatthew G. Knepley 87629d8e7bcSMatthew G. Knepley ierr = PetscHeaderCreate(s, _p_PetscSpace, struct _PetscSpaceOps, PETSCSPACE_CLASSID, "PetscSpace", "Linear Space", "PetscSpace", comm, PetscSpaceDestroy, PetscSpaceView);CHKERRQ(ierr); 87729d8e7bcSMatthew G. Knepley ierr = PetscMemzero(s->ops, sizeof(struct _PetscSpaceOps));CHKERRQ(ierr); 87829d8e7bcSMatthew G. Knepley 87929d8e7bcSMatthew G. Knepley s->order = 0; 88029d8e7bcSMatthew G. Knepley ierr = DMShellCreate(comm, &s->dm);CHKERRQ(ierr); 88129d8e7bcSMatthew G. Knepley 88229d8e7bcSMatthew G. Knepley *sp = s; 88329d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 88429d8e7bcSMatthew G. Knepley } 88529d8e7bcSMatthew G. Knepley 88629d8e7bcSMatthew G. Knepley #undef __FUNCT__ 88729d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceGetDimension" 88829d8e7bcSMatthew G. Knepley /* Dimension of the space, i.e. number of basis vectors */ 88929d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceGetDimension(PetscSpace sp, PetscInt *dim) 89029d8e7bcSMatthew G. Knepley { 89129d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 89229d8e7bcSMatthew G. Knepley 89329d8e7bcSMatthew G. Knepley PetscFunctionBegin; 89429d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 89529d8e7bcSMatthew G. Knepley PetscValidPointer(dim, 2); 89629d8e7bcSMatthew G. Knepley *dim = 0; 89729d8e7bcSMatthew G. Knepley if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 89829d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 89929d8e7bcSMatthew G. Knepley } 90029d8e7bcSMatthew G. Knepley 90129d8e7bcSMatthew G. Knepley #undef __FUNCT__ 90229d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceGetOrder" 90329d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceGetOrder(PetscSpace sp, PetscInt *order) 90429d8e7bcSMatthew G. Knepley { 90529d8e7bcSMatthew G. Knepley PetscFunctionBegin; 90629d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 90729d8e7bcSMatthew G. Knepley PetscValidPointer(order, 2); 90829d8e7bcSMatthew G. Knepley *order = sp->order; 90929d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 91029d8e7bcSMatthew G. Knepley } 91129d8e7bcSMatthew G. Knepley 91229d8e7bcSMatthew G. Knepley #undef __FUNCT__ 91329d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetOrder" 91429d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetOrder(PetscSpace sp, PetscInt order) 91529d8e7bcSMatthew G. Knepley { 91629d8e7bcSMatthew G. Knepley PetscFunctionBegin; 91729d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 91829d8e7bcSMatthew G. Knepley sp->order = order; 91929d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 92029d8e7bcSMatthew G. Knepley } 92129d8e7bcSMatthew G. Knepley 92229d8e7bcSMatthew G. Knepley #undef __FUNCT__ 92329d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceEvaluate" 92429d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceEvaluate(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 92529d8e7bcSMatthew G. Knepley { 92629d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 92729d8e7bcSMatthew G. Knepley 92829d8e7bcSMatthew G. Knepley PetscFunctionBegin; 92929d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 93029d8e7bcSMatthew G. Knepley PetscValidPointer(points, 3); 93129d8e7bcSMatthew G. Knepley if (B) PetscValidPointer(B, 4); 93229d8e7bcSMatthew G. Knepley if (D) PetscValidPointer(D, 5); 93329d8e7bcSMatthew G. Knepley if (H) PetscValidPointer(H, 6); 93429d8e7bcSMatthew G. Knepley if (sp->ops->evaluate) {ierr = (*sp->ops->evaluate)(sp, npoints, points, B, D, H);CHKERRQ(ierr);} 93529d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 93629d8e7bcSMatthew G. Knepley } 93729d8e7bcSMatthew G. Knepley 93829d8e7bcSMatthew G. Knepley #undef __FUNCT__ 93929d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetFromOptions_Polynomial" 94029d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscSpace sp) 94129d8e7bcSMatthew G. Knepley { 94229d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 94329d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 94429d8e7bcSMatthew G. Knepley 94529d8e7bcSMatthew G. Knepley PetscFunctionBegin; 94629d8e7bcSMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 94729d8e7bcSMatthew 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); 94829d8e7bcSMatthew G. Knepley ierr = PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);CHKERRQ(ierr); 94929d8e7bcSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 95029d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 95129d8e7bcSMatthew G. Knepley } 95229d8e7bcSMatthew G. Knepley 95329d8e7bcSMatthew G. Knepley #undef __FUNCT__ 95429d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialView_Ascii" 95529d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer viewer) 95629d8e7bcSMatthew G. Knepley { 95729d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 95829d8e7bcSMatthew G. Knepley PetscViewerFormat format; 95929d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 96029d8e7bcSMatthew G. Knepley 96129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 96229d8e7bcSMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 96329d8e7bcSMatthew G. Knepley if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 96429d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr); 96529d8e7bcSMatthew G. Knepley } else { 96629d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr); 96729d8e7bcSMatthew G. Knepley } 96829d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 96929d8e7bcSMatthew G. Knepley } 97029d8e7bcSMatthew G. Knepley 97129d8e7bcSMatthew G. Knepley #undef __FUNCT__ 97229d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceView_Polynomial" 97329d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer) 97429d8e7bcSMatthew G. Knepley { 97529d8e7bcSMatthew G. Knepley PetscBool iascii; 97629d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 97729d8e7bcSMatthew G. Knepley 97829d8e7bcSMatthew G. Knepley PetscFunctionBegin; 97929d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 98029d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 98129d8e7bcSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 98229d8e7bcSMatthew G. Knepley if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);} 98329d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 98429d8e7bcSMatthew G. Knepley } 98529d8e7bcSMatthew G. Knepley 98629d8e7bcSMatthew G. Knepley #undef __FUNCT__ 98729d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceSetUp_Polynomial" 98829d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp) 98929d8e7bcSMatthew G. Knepley { 99029d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 99129d8e7bcSMatthew G. Knepley PetscInt ndegree = sp->order+1; 99229d8e7bcSMatthew G. Knepley PetscInt dim = poly->numVariables; 99329d8e7bcSMatthew G. Knepley PetscInt pdim, deg; 99429d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 99529d8e7bcSMatthew G. Knepley 99629d8e7bcSMatthew G. Knepley PetscFunctionBegin; 99729d8e7bcSMatthew G. Knepley ierr = PetscMalloc(ndegree * sizeof(PetscInt), &poly->degrees);CHKERRQ(ierr); 99829d8e7bcSMatthew G. Knepley for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg; 99929d8e7bcSMatthew G. Knepley ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 100029d8e7bcSMatthew G. Knepley /* Create A (pdim x ndegree * dim) */ 100129d8e7bcSMatthew G. Knepley ierr = PetscMalloc(pdim*ndegree*dim * sizeof(PetscReal), &poly->A);CHKERRQ(ierr); 100229d8e7bcSMatthew G. Knepley ierr = PetscMemzero(poly->A, pdim*ndegree*dim * sizeof(PetscReal));CHKERRQ(ierr); 100329d8e7bcSMatthew G. Knepley /* Hardcode P_1: Here we need a way to iterate through the basis */ 100429d8e7bcSMatthew G. Knepley poly->A[(0*ndegree + 0)*dim + 0] = 1.0; /* 1 */ 100529d8e7bcSMatthew G. Knepley poly->A[(1*ndegree + 1)*dim + 0] = 1.0; /* x */ 100629d8e7bcSMatthew G. Knepley poly->A[(2*ndegree + 1)*dim + 1] = 1.0; /* y */ 100729d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 100829d8e7bcSMatthew G. Knepley } 100929d8e7bcSMatthew G. Knepley 101029d8e7bcSMatthew G. Knepley #undef __FUNCT__ 101129d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceDestroy_Polynomial" 101229d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp) 101329d8e7bcSMatthew G. Knepley { 101429d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 101529d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 101629d8e7bcSMatthew G. Knepley 101729d8e7bcSMatthew G. Knepley PetscFunctionBegin; 101829d8e7bcSMatthew G. Knepley ierr = PetscFree(poly->A);CHKERRQ(ierr); 101929d8e7bcSMatthew G. Knepley ierr = PetscFree(poly->degrees);CHKERRQ(ierr); 102029d8e7bcSMatthew G. Knepley ierr = PetscFree(poly);CHKERRQ(ierr); 102129d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 102229d8e7bcSMatthew G. Knepley } 102329d8e7bcSMatthew G. Knepley 102429d8e7bcSMatthew G. Knepley #undef __FUNCT__ 102529d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceGetDimension_Polynomial" 102629d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim) 102729d8e7bcSMatthew G. Knepley { 102829d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 102929d8e7bcSMatthew G. Knepley PetscInt deg = sp->order; 103029d8e7bcSMatthew G. Knepley PetscInt n = poly->numVariables, i; 103129d8e7bcSMatthew G. Knepley PetscReal D = 1.0; 103229d8e7bcSMatthew G. Knepley 103329d8e7bcSMatthew G. Knepley PetscFunctionBegin; 103429d8e7bcSMatthew G. Knepley for (i = 1; i <= n; ++i) { 103529d8e7bcSMatthew G. Knepley D *= ((PetscReal) (deg+i))/i; 103629d8e7bcSMatthew G. Knepley } 103729d8e7bcSMatthew G. Knepley *dim = (PetscInt) (D + 0.5); 103829d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 103929d8e7bcSMatthew G. Knepley } 104029d8e7bcSMatthew G. Knepley 104129d8e7bcSMatthew G. Knepley #undef __FUNCT__ 10423dc9baa2SMatthew G. Knepley #define __FUNCT__ "LatticePoint_Internal" 10433dc9baa2SMatthew G. Knepley /* 10443dc9baa2SMatthew G. Knepley LatticePoint_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to 'sum'. 10453dc9baa2SMatthew G. Knepley 10463dc9baa2SMatthew G. Knepley Input Parameters: 10473dc9baa2SMatthew G. Knepley + len - The length of the tuple 10483dc9baa2SMatthew G. Knepley . sum - The sum of all entries in the tuple 10493dc9baa2SMatthew G. Knepley - ind - The current multi-index of the tuple, initialized to the 0 tuple 10503dc9baa2SMatthew G. Knepley 10513dc9baa2SMatthew G. Knepley Output Parameter: 10523dc9baa2SMatthew G. Knepley + ind - The multi-index of the tuple, -1 indicates the iteration has terminated 10533dc9baa2SMatthew G. Knepley . tup - A tuple of len integers addig to sum 10543dc9baa2SMatthew G. Knepley 10553dc9baa2SMatthew G. Knepley Level: developer 10563dc9baa2SMatthew G. Knepley 10573dc9baa2SMatthew G. Knepley .seealso: 10583dc9baa2SMatthew G. Knepley */ 10593dc9baa2SMatthew G. Knepley static PetscErrorCode LatticePoint_Internal(PetscInt len, PetscInt sum, PetscInt ind[], PetscInt tup[]) 10603dc9baa2SMatthew G. Knepley { 10613dc9baa2SMatthew G. Knepley PetscInt i; 10623dc9baa2SMatthew G. Knepley PetscErrorCode ierr; 10633dc9baa2SMatthew G. Knepley 10643dc9baa2SMatthew G. Knepley PetscFunctionBegin; 10653dc9baa2SMatthew G. Knepley if (len == 1) { 10663dc9baa2SMatthew G. Knepley ind[0] = -1; 10673dc9baa2SMatthew G. Knepley tup[0] = sum; 10683dc9baa2SMatthew G. Knepley } else if (sum == 0) { 10693dc9baa2SMatthew G. Knepley for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;} 10703dc9baa2SMatthew G. Knepley } else { 10713dc9baa2SMatthew G. Knepley tup[0] = sum - ind[0]; 10723dc9baa2SMatthew G. Knepley ierr = LatticePoint_Internal(len-1, ind[0], &ind[1], &tup[1]);CHKERRQ(ierr); 10733dc9baa2SMatthew G. Knepley if (ind[1] < 0) { 10743dc9baa2SMatthew G. Knepley if (ind[0] == sum) {ind[0] = -1;} 10753dc9baa2SMatthew G. Knepley else {ind[1] = 0; ++ind[0];} 10763dc9baa2SMatthew G. Knepley } 10773dc9baa2SMatthew G. Knepley } 10783dc9baa2SMatthew G. Knepley PetscFunctionReturn(0); 10793dc9baa2SMatthew G. Knepley } 10803dc9baa2SMatthew G. Knepley 10813dc9baa2SMatthew G. Knepley #undef __FUNCT__ 108229d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceEvaluate_Polynomial" 108329d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 108429d8e7bcSMatthew G. Knepley { 108529d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 108629d8e7bcSMatthew G. Knepley DM dm = sp->dm; 108729d8e7bcSMatthew G. Knepley PetscInt ndegree = sp->order+1; 108829d8e7bcSMatthew G. Knepley PetscInt *degrees = poly->degrees; 108929d8e7bcSMatthew G. Knepley PetscInt dim = poly->numVariables; 109029d8e7bcSMatthew G. Knepley PetscReal *A = poly->A; 109129d8e7bcSMatthew G. Knepley PetscReal *lpoints, *tmp, *LB, *LD, *LH; 10923dc9baa2SMatthew G. Knepley PetscInt *ind, *tup; 10933dc9baa2SMatthew G. Knepley PetscInt pdim, d, i, p, deg, o; 109429d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 109529d8e7bcSMatthew G. Knepley 109629d8e7bcSMatthew G. Knepley PetscFunctionBegin; 109729d8e7bcSMatthew G. Knepley ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 109829d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr); 109929d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr); 110029d8e7bcSMatthew G. Knepley if (B) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);} 110129d8e7bcSMatthew G. Knepley if (D) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);} 110229d8e7bcSMatthew G. Knepley if (H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);} 110329d8e7bcSMatthew G. Knepley for (d = 0; d < dim; ++d) { 110429d8e7bcSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 110529d8e7bcSMatthew G. Knepley lpoints[p] = points[p*dim+d]; 110629d8e7bcSMatthew G. Knepley } 110729d8e7bcSMatthew G. Knepley ierr = PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);CHKERRQ(ierr); 110829d8e7bcSMatthew G. Knepley /* LB, LD, LH (ndegree * dim x npoints) */ 110929d8e7bcSMatthew G. Knepley for (deg = 0; deg < ndegree; ++deg) { 111029d8e7bcSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 111129d8e7bcSMatthew G. Knepley if (B) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg]; 111229d8e7bcSMatthew G. Knepley if (D) LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg]; 111329d8e7bcSMatthew G. Knepley if (H) LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg]; 111429d8e7bcSMatthew G. Knepley } 111529d8e7bcSMatthew G. Knepley } 111629d8e7bcSMatthew G. Knepley } 111729d8e7bcSMatthew G. Knepley /* Multiply by A (pdim x ndegree * dim) */ 11183dc9baa2SMatthew G. Knepley ierr = PetscMalloc2(dim,PetscInt,&ind,dim,PetscInt,&tup);CHKERRQ(ierr); 111929d8e7bcSMatthew G. Knepley if (B) { 11203dc9baa2SMatthew G. Knepley /* B (npoints x pdim) */ 11213dc9baa2SMatthew G. Knepley i = 0; 11223dc9baa2SMatthew G. Knepley for (o = 0; o <= sp->order; ++o) { 11233dc9baa2SMatthew G. Knepley ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr); 11243dc9baa2SMatthew G. Knepley while (ind[0] >= 0) { 11253dc9baa2SMatthew G. Knepley ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr); 112629d8e7bcSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 11273dc9baa2SMatthew G. Knepley B[p*pdim + i] = 1.0; 112829d8e7bcSMatthew G. Knepley for (d = 0; d < dim; ++d) { 11293dc9baa2SMatthew G. Knepley B[p*pdim + i] *= LB[(tup[d]*dim + d)*npoints + p]; 11303dc9baa2SMatthew G. Knepley } 11313dc9baa2SMatthew G. Knepley } 11323dc9baa2SMatthew G. Knepley ++i; 113329d8e7bcSMatthew G. Knepley } 113429d8e7bcSMatthew G. Knepley } 113529d8e7bcSMatthew G. Knepley } 11363dc9baa2SMatthew G. Knepley ierr = PetscFree2(ind,tup);CHKERRQ(ierr); 113729d8e7bcSMatthew G. Knepley if (D) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code first derivatives"); 113829d8e7bcSMatthew G. Knepley if (H) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code second derivatives"); 113929d8e7bcSMatthew G. Knepley if (B) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);} 114029d8e7bcSMatthew G. Knepley if (D) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);} 114129d8e7bcSMatthew G. Knepley if (H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);} 114229d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr); 114329d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr); 114429d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 114529d8e7bcSMatthew G. Knepley } 114629d8e7bcSMatthew G. Knepley 114729d8e7bcSMatthew G. Knepley #undef __FUNCT__ 114829d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceInitialize_Polynomial" 114929d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp) 115029d8e7bcSMatthew G. Knepley { 115129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 115229d8e7bcSMatthew G. Knepley sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial; 115329d8e7bcSMatthew G. Knepley sp->ops->setup = PetscSpaceSetUp_Polynomial; 115429d8e7bcSMatthew G. Knepley sp->ops->view = PetscSpaceView_Polynomial; 115529d8e7bcSMatthew G. Knepley sp->ops->destroy = PetscSpaceDestroy_Polynomial; 115629d8e7bcSMatthew G. Knepley sp->ops->getdimension = PetscSpaceGetDimension_Polynomial; 115729d8e7bcSMatthew G. Knepley sp->ops->evaluate = PetscSpaceEvaluate_Polynomial; 115829d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 115929d8e7bcSMatthew G. Knepley } 116029d8e7bcSMatthew G. Knepley 116129d8e7bcSMatthew G. Knepley /*MC 116229d8e7bcSMatthew G. Knepley PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of linear polynomials. 116329d8e7bcSMatthew G. Knepley 116429d8e7bcSMatthew G. Knepley Level: intermediate 116529d8e7bcSMatthew G. Knepley 116629d8e7bcSMatthew G. Knepley .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType() 116729d8e7bcSMatthew G. Knepley M*/ 116829d8e7bcSMatthew G. Knepley 116929d8e7bcSMatthew G. Knepley #undef __FUNCT__ 117029d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpaceCreate_Polynomial" 117129d8e7bcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp) 117229d8e7bcSMatthew G. Knepley { 117329d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly; 117429d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 117529d8e7bcSMatthew G. Knepley 117629d8e7bcSMatthew G. Knepley PetscFunctionBegin; 117729d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 117829d8e7bcSMatthew G. Knepley ierr = PetscNewLog(sp, PetscSpace_Poly, &poly);CHKERRQ(ierr); 117929d8e7bcSMatthew G. Knepley sp->data = poly; 118029d8e7bcSMatthew G. Knepley 118129d8e7bcSMatthew G. Knepley poly->numVariables = 0; 118229d8e7bcSMatthew G. Knepley poly->symmetric = PETSC_FALSE; 118329d8e7bcSMatthew G. Knepley poly->degrees = NULL; 118429d8e7bcSMatthew G. Knepley poly->A = NULL; 118529d8e7bcSMatthew G. Knepley 118629d8e7bcSMatthew G. Knepley ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr); 118729d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 118829d8e7bcSMatthew G. Knepley } 118929d8e7bcSMatthew G. Knepley 119029d8e7bcSMatthew G. Knepley #undef __FUNCT__ 119129d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialSetSymmetric" 119229d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool sym) 119329d8e7bcSMatthew G. Knepley { 119429d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 119529d8e7bcSMatthew G. Knepley 119629d8e7bcSMatthew G. Knepley PetscFunctionBegin; 119729d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 119829d8e7bcSMatthew G. Knepley poly->symmetric = sym; 119929d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 120029d8e7bcSMatthew G. Knepley } 120129d8e7bcSMatthew G. Knepley 120229d8e7bcSMatthew G. Knepley #undef __FUNCT__ 120329d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialGetSymmetric" 120429d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace sp, PetscBool *sym) 120529d8e7bcSMatthew G. Knepley { 120629d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 120729d8e7bcSMatthew G. Knepley 120829d8e7bcSMatthew G. Knepley PetscFunctionBegin; 120929d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 121029d8e7bcSMatthew G. Knepley PetscValidPointer(sym, 2); 121129d8e7bcSMatthew G. Knepley *sym = poly->symmetric; 121229d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 121329d8e7bcSMatthew G. Knepley } 121429d8e7bcSMatthew G. Knepley 121529d8e7bcSMatthew G. Knepley #undef __FUNCT__ 121629d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialSetNumVariables" 121729d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace sp, PetscInt n) 121829d8e7bcSMatthew G. Knepley { 121929d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 122029d8e7bcSMatthew G. Knepley 122129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 122229d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 122329d8e7bcSMatthew G. Knepley poly->numVariables = n; 122429d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 122529d8e7bcSMatthew G. Knepley } 122629d8e7bcSMatthew G. Knepley 122729d8e7bcSMatthew G. Knepley #undef __FUNCT__ 122829d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscSpacePolynomialGetNumVariables" 122929d8e7bcSMatthew G. Knepley PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace sp, PetscInt *n) 123029d8e7bcSMatthew G. Knepley { 123129d8e7bcSMatthew G. Knepley PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 123229d8e7bcSMatthew G. Knepley 123329d8e7bcSMatthew G. Knepley PetscFunctionBegin; 123429d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 123529d8e7bcSMatthew G. Knepley PetscValidPointer(n, 2); 123629d8e7bcSMatthew G. Knepley *n = poly->numVariables; 123729d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 123829d8e7bcSMatthew G. Knepley } 123929d8e7bcSMatthew G. Knepley 124029d8e7bcSMatthew G. Knepley 124129d8e7bcSMatthew G. Knepley PetscInt PETSCDUALSPACE_CLASSID = 0; 124229d8e7bcSMatthew G. Knepley 124329d8e7bcSMatthew G. Knepley PetscFunctionList PetscDualSpaceList = NULL; 124429d8e7bcSMatthew G. Knepley PetscBool PetscDualSpaceRegisterAllCalled = PETSC_FALSE; 124529d8e7bcSMatthew G. Knepley 124629d8e7bcSMatthew G. Knepley #undef __FUNCT__ 124729d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceRegister" 124829d8e7bcSMatthew G. Knepley /*@C 124929d8e7bcSMatthew G. Knepley PetscDualSpaceRegister - Adds a new PetscDualSpace implementation 125029d8e7bcSMatthew G. Knepley 125129d8e7bcSMatthew G. Knepley Not Collective 125229d8e7bcSMatthew G. Knepley 125329d8e7bcSMatthew G. Knepley Input Parameters: 125429d8e7bcSMatthew G. Knepley + name - The name of a new user-defined creation routine 125529d8e7bcSMatthew G. Knepley - create_func - The creation routine itself 125629d8e7bcSMatthew G. Knepley 125729d8e7bcSMatthew G. Knepley Notes: 125829d8e7bcSMatthew G. Knepley PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces 125929d8e7bcSMatthew G. Knepley 126029d8e7bcSMatthew G. Knepley Sample usage: 126129d8e7bcSMatthew G. Knepley .vb 126229d8e7bcSMatthew G. Knepley PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate); 126329d8e7bcSMatthew G. Knepley .ve 126429d8e7bcSMatthew G. Knepley 126529d8e7bcSMatthew G. Knepley Then, your PetscDualSpace type can be chosen with the procedural interface via 126629d8e7bcSMatthew G. Knepley .vb 126729d8e7bcSMatthew G. Knepley PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 126829d8e7bcSMatthew G. Knepley PetscDualSpaceSetType(PetscDualSpace, "my_dual_space"); 126929d8e7bcSMatthew G. Knepley .ve 127029d8e7bcSMatthew G. Knepley or at runtime via the option 127129d8e7bcSMatthew G. Knepley .vb 127229d8e7bcSMatthew G. Knepley -petscdualspace_type my_dual_space 127329d8e7bcSMatthew G. Knepley .ve 127429d8e7bcSMatthew G. Knepley 127529d8e7bcSMatthew G. Knepley Level: advanced 127629d8e7bcSMatthew G. Knepley 127729d8e7bcSMatthew G. Knepley .keywords: PetscDualSpace, register 127829d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy() 127929d8e7bcSMatthew G. Knepley 128029d8e7bcSMatthew G. Knepley @*/ 128129d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace)) 128229d8e7bcSMatthew G. Knepley { 128329d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 128429d8e7bcSMatthew G. Knepley 128529d8e7bcSMatthew G. Knepley PetscFunctionBegin; 128629d8e7bcSMatthew G. Knepley ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr); 128729d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 128829d8e7bcSMatthew G. Knepley } 128929d8e7bcSMatthew G. Knepley 129029d8e7bcSMatthew G. Knepley #undef __FUNCT__ 129129d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetType" 129229d8e7bcSMatthew G. Knepley /*@C 129329d8e7bcSMatthew G. Knepley PetscDualSpaceSetType - Builds a particular PetscDualSpace 129429d8e7bcSMatthew G. Knepley 129529d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 129629d8e7bcSMatthew G. Knepley 129729d8e7bcSMatthew G. Knepley Input Parameters: 129829d8e7bcSMatthew G. Knepley + sp - The PetscDualSpace object 129929d8e7bcSMatthew G. Knepley - name - The kind of space 130029d8e7bcSMatthew G. Knepley 130129d8e7bcSMatthew G. Knepley Options Database Key: 130229d8e7bcSMatthew G. Knepley . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types 130329d8e7bcSMatthew G. Knepley 130429d8e7bcSMatthew G. Knepley Level: intermediate 130529d8e7bcSMatthew G. Knepley 130629d8e7bcSMatthew G. Knepley .keywords: PetscDualSpace, set, type 130729d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate() 130829d8e7bcSMatthew G. Knepley @*/ 130929d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name) 131029d8e7bcSMatthew G. Knepley { 131129d8e7bcSMatthew G. Knepley PetscErrorCode (*r)(PetscDualSpace); 131229d8e7bcSMatthew G. Knepley PetscBool match; 131329d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 131429d8e7bcSMatthew G. Knepley 131529d8e7bcSMatthew G. Knepley PetscFunctionBegin; 131629d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 131729d8e7bcSMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 131829d8e7bcSMatthew G. Knepley if (match) PetscFunctionReturn(0); 131929d8e7bcSMatthew G. Knepley 132029d8e7bcSMatthew G. Knepley if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);} 132129d8e7bcSMatthew G. Knepley ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr); 132229d8e7bcSMatthew G. Knepley if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name); 132329d8e7bcSMatthew G. Knepley 132429d8e7bcSMatthew G. Knepley if (sp->ops->destroy) { 132529d8e7bcSMatthew G. Knepley ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 132629d8e7bcSMatthew G. Knepley sp->ops->destroy = NULL; 132729d8e7bcSMatthew G. Knepley } 132829d8e7bcSMatthew G. Knepley ierr = (*r)(sp);CHKERRQ(ierr); 132929d8e7bcSMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 133029d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 133129d8e7bcSMatthew G. Knepley } 133229d8e7bcSMatthew G. Knepley 133329d8e7bcSMatthew G. Knepley #undef __FUNCT__ 133429d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetType" 133529d8e7bcSMatthew G. Knepley /*@C 133629d8e7bcSMatthew G. Knepley PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object. 133729d8e7bcSMatthew G. Knepley 133829d8e7bcSMatthew G. Knepley Not Collective 133929d8e7bcSMatthew G. Knepley 134029d8e7bcSMatthew G. Knepley Input Parameter: 134129d8e7bcSMatthew G. Knepley . dm - The PetscDualSpace 134229d8e7bcSMatthew G. Knepley 134329d8e7bcSMatthew G. Knepley Output Parameter: 134429d8e7bcSMatthew G. Knepley . name - The PetscDualSpace type name 134529d8e7bcSMatthew G. Knepley 134629d8e7bcSMatthew G. Knepley Level: intermediate 134729d8e7bcSMatthew G. Knepley 134829d8e7bcSMatthew G. Knepley .keywords: PetscDualSpace, get, type, name 134929d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate() 135029d8e7bcSMatthew G. Knepley @*/ 135129d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name) 135229d8e7bcSMatthew G. Knepley { 135329d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 135429d8e7bcSMatthew G. Knepley 135529d8e7bcSMatthew G. Knepley PetscFunctionBegin; 135629d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 135729d8e7bcSMatthew G. Knepley PetscValidCharPointer(name, 2); 135829d8e7bcSMatthew G. Knepley if (!PetscDualSpaceRegisterAllCalled) { 135929d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr); 136029d8e7bcSMatthew G. Knepley } 136129d8e7bcSMatthew G. Knepley *name = ((PetscObject) sp)->type_name; 136229d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 136329d8e7bcSMatthew G. Knepley } 136429d8e7bcSMatthew G. Knepley 136529d8e7bcSMatthew G. Knepley #undef __FUNCT__ 136629d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceView" 136729d8e7bcSMatthew G. Knepley /*@C 136829d8e7bcSMatthew G. Knepley PetscDualSpaceView - Views a PetscDualSpace 136929d8e7bcSMatthew G. Knepley 137029d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 137129d8e7bcSMatthew G. Knepley 137229d8e7bcSMatthew G. Knepley Input Parameter: 137329d8e7bcSMatthew G. Knepley + sp - the PetscDualSpace object to view 137429d8e7bcSMatthew G. Knepley - v - the viewer 137529d8e7bcSMatthew G. Knepley 137629d8e7bcSMatthew G. Knepley Level: developer 137729d8e7bcSMatthew G. Knepley 137829d8e7bcSMatthew G. Knepley .seealso PetscDualSpaceDestroy() 137929d8e7bcSMatthew G. Knepley @*/ 138029d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v) 138129d8e7bcSMatthew G. Knepley { 138229d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 138329d8e7bcSMatthew G. Knepley 138429d8e7bcSMatthew G. Knepley PetscFunctionBegin; 138529d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 138629d8e7bcSMatthew G. Knepley if (!v) { 138729d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr); 138829d8e7bcSMatthew G. Knepley } 138929d8e7bcSMatthew G. Knepley if (sp->ops->view) { 139029d8e7bcSMatthew G. Knepley ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr); 139129d8e7bcSMatthew G. Knepley } 139229d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 139329d8e7bcSMatthew G. Knepley } 139429d8e7bcSMatthew G. Knepley 139529d8e7bcSMatthew G. Knepley #undef __FUNCT__ 139629d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceViewFromOptions" 139729d8e7bcSMatthew G. Knepley /* 139829d8e7bcSMatthew G. Knepley PetscDualSpaceViewFromOptions - Processes command line options to determine if/how a PetscDualSpace is to be viewed. 139929d8e7bcSMatthew G. Knepley 140029d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 140129d8e7bcSMatthew G. Knepley 140229d8e7bcSMatthew G. Knepley Input Parameters: 140329d8e7bcSMatthew G. Knepley + sp - the PetscDualSpace 140429d8e7bcSMatthew G. Knepley . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd' 140529d8e7bcSMatthew G. Knepley - optionname - option to activate viewing 140629d8e7bcSMatthew G. Knepley 140729d8e7bcSMatthew G. Knepley Level: intermediate 140829d8e7bcSMatthew G. Knepley 140929d8e7bcSMatthew G. Knepley .keywords: PetscDualSpace, view, options, database 141029d8e7bcSMatthew G. Knepley .seealso: VecViewFromOptions(), MatViewFromOptions() 141129d8e7bcSMatthew G. Knepley */ 141229d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace sp, const char prefix[], const char optionname[]) 141329d8e7bcSMatthew G. Knepley { 141429d8e7bcSMatthew G. Knepley PetscViewer viewer; 141529d8e7bcSMatthew G. Knepley PetscViewerFormat format; 141629d8e7bcSMatthew G. Knepley PetscBool flg; 141729d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 141829d8e7bcSMatthew G. Knepley 141929d8e7bcSMatthew G. Knepley PetscFunctionBegin; 142029d8e7bcSMatthew G. Knepley if (prefix) { 142129d8e7bcSMatthew G. Knepley ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 142229d8e7bcSMatthew G. Knepley } else { 142329d8e7bcSMatthew G. Knepley ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr); 142429d8e7bcSMatthew G. Knepley } 142529d8e7bcSMatthew G. Knepley if (flg) { 142629d8e7bcSMatthew G. Knepley ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 142729d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceView(sp, viewer);CHKERRQ(ierr); 142829d8e7bcSMatthew G. Knepley ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 142929d8e7bcSMatthew G. Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 143029d8e7bcSMatthew G. Knepley } 143129d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 143229d8e7bcSMatthew G. Knepley } 143329d8e7bcSMatthew G. Knepley 143429d8e7bcSMatthew G. Knepley #undef __FUNCT__ 143529d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetFromOptions" 143629d8e7bcSMatthew G. Knepley /*@ 143729d8e7bcSMatthew G. Knepley PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database 143829d8e7bcSMatthew G. Knepley 143929d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 144029d8e7bcSMatthew G. Knepley 144129d8e7bcSMatthew G. Knepley Input Parameter: 144229d8e7bcSMatthew G. Knepley . sp - the PetscDualSpace object to set options for 144329d8e7bcSMatthew G. Knepley 144429d8e7bcSMatthew G. Knepley Options Database: 144529d8e7bcSMatthew G. Knepley . -petscspace_order the approximation order of the space 144629d8e7bcSMatthew G. Knepley 144729d8e7bcSMatthew G. Knepley Level: developer 144829d8e7bcSMatthew G. Knepley 144929d8e7bcSMatthew G. Knepley .seealso PetscDualSpaceView() 145029d8e7bcSMatthew G. Knepley @*/ 145129d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp) 145229d8e7bcSMatthew G. Knepley { 145329d8e7bcSMatthew G. Knepley const char *defaultType; 145429d8e7bcSMatthew G. Knepley char typename[256]; 145529d8e7bcSMatthew G. Knepley PetscBool flg; 145629d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 145729d8e7bcSMatthew G. Knepley 145829d8e7bcSMatthew G. Knepley PetscFunctionBegin; 145929d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 146029d8e7bcSMatthew G. Knepley if (!((PetscObject) sp)->type_name) { 146129d8e7bcSMatthew G. Knepley defaultType = PETSCDUALSPACELAGRANGE; 146229d8e7bcSMatthew G. Knepley } else { 146329d8e7bcSMatthew G. Knepley defaultType = ((PetscObject) sp)->type_name; 146429d8e7bcSMatthew G. Knepley } 146529d8e7bcSMatthew G. Knepley if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 146629d8e7bcSMatthew G. Knepley 146729d8e7bcSMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 146829d8e7bcSMatthew G. Knepley ierr = PetscOptionsList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr); 146929d8e7bcSMatthew G. Knepley if (flg) { 147029d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceSetType(sp, typename);CHKERRQ(ierr); 147129d8e7bcSMatthew G. Knepley } else if (!((PetscObject) sp)->type_name) { 147229d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr); 147329d8e7bcSMatthew G. Knepley } 147429d8e7bcSMatthew G. Knepley ierr = PetscOptionsInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr); 147529d8e7bcSMatthew G. Knepley if (sp->ops->setfromoptions) { 147629d8e7bcSMatthew G. Knepley ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr); 147729d8e7bcSMatthew G. Knepley } 147829d8e7bcSMatthew G. Knepley /* process any options handlers added with PetscObjectAddOptionsHandler() */ 147929d8e7bcSMatthew G. Knepley ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr); 148029d8e7bcSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 148129d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr); 148229d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 148329d8e7bcSMatthew G. Knepley } 148429d8e7bcSMatthew G. Knepley 148529d8e7bcSMatthew G. Knepley #undef __FUNCT__ 148629d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetUp" 148729d8e7bcSMatthew G. Knepley /*@C 148829d8e7bcSMatthew G. Knepley PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace 148929d8e7bcSMatthew G. Knepley 149029d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 149129d8e7bcSMatthew G. Knepley 149229d8e7bcSMatthew G. Knepley Input Parameter: 149329d8e7bcSMatthew G. Knepley . sp - the PetscDualSpace object to setup 149429d8e7bcSMatthew G. Knepley 149529d8e7bcSMatthew G. Knepley Level: developer 149629d8e7bcSMatthew G. Knepley 149729d8e7bcSMatthew G. Knepley .seealso PetscDualSpaceView(), PetscDualSpaceDestroy() 149829d8e7bcSMatthew G. Knepley @*/ 149929d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp) 150029d8e7bcSMatthew G. Knepley { 150129d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 150229d8e7bcSMatthew G. Knepley 150329d8e7bcSMatthew G. Knepley PetscFunctionBegin; 150429d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 150529d8e7bcSMatthew G. Knepley if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 150629d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 150729d8e7bcSMatthew G. Knepley } 150829d8e7bcSMatthew G. Knepley 150929d8e7bcSMatthew G. Knepley #undef __FUNCT__ 151029d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceDestroy" 151129d8e7bcSMatthew G. Knepley /*@ 151229d8e7bcSMatthew G. Knepley PetscDualSpaceDestroy - Destroys a PetscDualSpace object 151329d8e7bcSMatthew G. Knepley 151429d8e7bcSMatthew G. Knepley Collective on PetscDualSpace 151529d8e7bcSMatthew G. Knepley 151629d8e7bcSMatthew G. Knepley Input Parameter: 151729d8e7bcSMatthew G. Knepley . sp - the PetscDualSpace object to destroy 151829d8e7bcSMatthew G. Knepley 151929d8e7bcSMatthew G. Knepley Level: developer 152029d8e7bcSMatthew G. Knepley 152129d8e7bcSMatthew G. Knepley .seealso PetscDualSpaceView() 152229d8e7bcSMatthew G. Knepley @*/ 152329d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp) 152429d8e7bcSMatthew G. Knepley { 152529d8e7bcSMatthew G. Knepley PetscInt dim, f; 152629d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 152729d8e7bcSMatthew G. Knepley 152829d8e7bcSMatthew G. Knepley PetscFunctionBegin; 152929d8e7bcSMatthew G. Knepley if (!*sp) PetscFunctionReturn(0); 153029d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1); 153129d8e7bcSMatthew G. Knepley 153229d8e7bcSMatthew G. Knepley if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 153329d8e7bcSMatthew G. Knepley ((PetscObject) (*sp))->refct = 0; 153429d8e7bcSMatthew G. Knepley /* if memory was published with AMS then destroy it */ 153529d8e7bcSMatthew G. Knepley ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr); 153629d8e7bcSMatthew G. Knepley 153754da605fSMatthew G. Knepley ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr); 153829d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr); 153929d8e7bcSMatthew G. Knepley for (f = 0; f < dim; ++f) { 154029d8e7bcSMatthew G. Knepley /* ierr = PetscQuadratureDestroy((*sp)->functional[f]);CHKERRQ(ierr); */ 154129d8e7bcSMatthew G. Knepley } 154229d8e7bcSMatthew G. Knepley ierr = PetscFree((*sp)->functional);CHKERRQ(ierr); 154329d8e7bcSMatthew G. Knepley 154429d8e7bcSMatthew G. Knepley ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr); 154529d8e7bcSMatthew G. Knepley ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 154629d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 154729d8e7bcSMatthew G. Knepley } 154829d8e7bcSMatthew G. Knepley 154929d8e7bcSMatthew G. Knepley #undef __FUNCT__ 155029d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceCreate" 155129d8e7bcSMatthew G. Knepley /*@ 155229d8e7bcSMatthew G. Knepley PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType(). 155329d8e7bcSMatthew G. Knepley 155429d8e7bcSMatthew G. Knepley Collective on MPI_Comm 155529d8e7bcSMatthew G. Knepley 155629d8e7bcSMatthew G. Knepley Input Parameter: 155729d8e7bcSMatthew G. Knepley . comm - The communicator for the PetscDualSpace object 155829d8e7bcSMatthew G. Knepley 155929d8e7bcSMatthew G. Knepley Output Parameter: 156029d8e7bcSMatthew G. Knepley . sp - The PetscDualSpace object 156129d8e7bcSMatthew G. Knepley 156229d8e7bcSMatthew G. Knepley Level: beginner 156329d8e7bcSMatthew G. Knepley 156429d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE 156529d8e7bcSMatthew G. Knepley @*/ 156629d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp) 156729d8e7bcSMatthew G. Knepley { 156829d8e7bcSMatthew G. Knepley PetscDualSpace s; 156929d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 157029d8e7bcSMatthew G. Knepley 157129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 157229d8e7bcSMatthew G. Knepley PetscValidPointer(sp, 2); 157329d8e7bcSMatthew G. Knepley *sp = NULL; 157429d8e7bcSMatthew G. Knepley #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 157529d8e7bcSMatthew G. Knepley ierr = PetscFEInitializePackage();CHKERRQ(ierr); 157629d8e7bcSMatthew G. Knepley #endif 157729d8e7bcSMatthew G. Knepley 157829d8e7bcSMatthew G. Knepley ierr = PetscHeaderCreate(s, _p_PetscDualSpace, struct _PetscDualSpaceOps, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr); 157929d8e7bcSMatthew G. Knepley ierr = PetscMemzero(s->ops, sizeof(struct _PetscDualSpaceOps));CHKERRQ(ierr); 158029d8e7bcSMatthew G. Knepley 158129d8e7bcSMatthew G. Knepley s->order = 0; 158229d8e7bcSMatthew G. Knepley 158329d8e7bcSMatthew G. Knepley *sp = s; 158429d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 158529d8e7bcSMatthew G. Knepley } 158629d8e7bcSMatthew G. Knepley 158729d8e7bcSMatthew G. Knepley #undef __FUNCT__ 158829d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetDM" 158929d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm) 159029d8e7bcSMatthew G. Knepley { 159129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 159229d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 159329d8e7bcSMatthew G. Knepley PetscValidPointer(dm, 2); 159429d8e7bcSMatthew G. Knepley *dm = sp->dm; 159529d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 159629d8e7bcSMatthew G. Knepley } 159729d8e7bcSMatthew G. Knepley 159829d8e7bcSMatthew G. Knepley #undef __FUNCT__ 159929d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetDM" 160029d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm) 160129d8e7bcSMatthew G. Knepley { 160254da605fSMatthew G. Knepley PetscErrorCode ierr; 160354da605fSMatthew G. Knepley 160429d8e7bcSMatthew G. Knepley PetscFunctionBegin; 160529d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 160629d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 160754da605fSMatthew G. Knepley ierr = DMDestroy(&sp->dm);CHKERRQ(ierr); 160854da605fSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 160929d8e7bcSMatthew G. Knepley sp->dm = dm; 161029d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 161129d8e7bcSMatthew G. Knepley } 161229d8e7bcSMatthew G. Knepley 161329d8e7bcSMatthew G. Knepley #undef __FUNCT__ 161429d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetOrder" 161529d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order) 161629d8e7bcSMatthew G. Knepley { 161729d8e7bcSMatthew G. Knepley PetscFunctionBegin; 161829d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 161929d8e7bcSMatthew G. Knepley PetscValidPointer(order, 2); 162029d8e7bcSMatthew G. Knepley *order = sp->order; 162129d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 162229d8e7bcSMatthew G. Knepley } 162329d8e7bcSMatthew G. Knepley 162429d8e7bcSMatthew G. Knepley #undef __FUNCT__ 162529d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetOrder" 162629d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order) 162729d8e7bcSMatthew G. Knepley { 162829d8e7bcSMatthew G. Knepley PetscFunctionBegin; 162929d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 163029d8e7bcSMatthew G. Knepley sp->order = order; 163129d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 163229d8e7bcSMatthew G. Knepley } 163329d8e7bcSMatthew G. Knepley 163429d8e7bcSMatthew G. Knepley #undef __FUNCT__ 163529d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetFunctional" 163629d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional) 163729d8e7bcSMatthew G. Knepley { 163829d8e7bcSMatthew G. Knepley PetscInt dim; 163929d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 164029d8e7bcSMatthew G. Knepley 164129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 164229d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 164329d8e7bcSMatthew G. Knepley PetscValidPointer(functional, 3); 164429d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 164529d8e7bcSMatthew 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); 164629d8e7bcSMatthew G. Knepley *functional = sp->functional[i]; 164729d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 164829d8e7bcSMatthew G. Knepley } 164929d8e7bcSMatthew G. Knepley 165029d8e7bcSMatthew G. Knepley #undef __FUNCT__ 165129d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetDimension" 165229d8e7bcSMatthew G. Knepley /* Dimension of the space, i.e. number of basis vectors */ 165329d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim) 165429d8e7bcSMatthew G. Knepley { 165529d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 165629d8e7bcSMatthew G. Knepley 165729d8e7bcSMatthew G. Knepley PetscFunctionBegin; 165829d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 165929d8e7bcSMatthew G. Knepley PetscValidPointer(dim, 2); 166029d8e7bcSMatthew G. Knepley *dim = 0; 166129d8e7bcSMatthew G. Knepley if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 166229d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 166329d8e7bcSMatthew G. Knepley } 166429d8e7bcSMatthew G. Knepley 166529d8e7bcSMatthew G. Knepley #undef __FUNCT__ 166629d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceSetUp_Lagrange" 166729d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 166829d8e7bcSMatthew G. Knepley { 166929d8e7bcSMatthew G. Knepley DM dm = sp->dm; 167029d8e7bcSMatthew G. Knepley PetscInt order = sp->order; 167129d8e7bcSMatthew G. Knepley PetscSection csection; 167229d8e7bcSMatthew G. Knepley Vec coordinates; 167359804f93SMatthew G. Knepley PetscReal *qpoints, *qweights; 167459804f93SMatthew G. Knepley PetscInt depth, dim, pdim, *pStart, *pEnd, coneSize, d, n, f = 0; 167529d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 167629d8e7bcSMatthew G. Knepley 167729d8e7bcSMatthew G. Knepley PetscFunctionBegin; 167829d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 167929d8e7bcSMatthew G. Knepley ierr = PetscMalloc(pdim * sizeof(PetscQuadrature), &sp->functional);CHKERRQ(ierr); 168029d8e7bcSMatthew G. Knepley /* Classify element type */ 168129d8e7bcSMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 168259804f93SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 168359804f93SMatthew G. Knepley ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr); 168459804f93SMatthew G. Knepley for (d = 0; d < depth; ++d) { 168559804f93SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 168659804f93SMatthew G. Knepley } 168759804f93SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, pStart[depth], &coneSize);CHKERRQ(ierr); 168829d8e7bcSMatthew G. Knepley ierr = DMPlexGetCoordinateSection(dm, &csection);CHKERRQ(ierr); 168929d8e7bcSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 169029d8e7bcSMatthew G. Knepley if (coneSize == dim+1) { 169129d8e7bcSMatthew G. Knepley PetscInt *closure = NULL, closureSize, c; 169229d8e7bcSMatthew G. Knepley 169329d8e7bcSMatthew G. Knepley /* Simplex */ 169459804f93SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, pStart[depth], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 169529d8e7bcSMatthew G. Knepley for (c = 0; c < closureSize*2; c += 2) { 169629d8e7bcSMatthew G. Knepley const PetscInt p = closure[c]; 169759804f93SMatthew G. Knepley 169859804f93SMatthew G. Knepley if ((p >= pStart[0]) && (p < pEnd[0])) { 169929d8e7bcSMatthew G. Knepley /* Vertices */ 1700*d9dbb98eSMatthew G. Knepley const PetscScalar *coords; 170129d8e7bcSMatthew G. Knepley PetscInt dof, off, d; 170229d8e7bcSMatthew G. Knepley 170359804f93SMatthew G. Knepley if (order < 1) continue; 170429d8e7bcSMatthew G. Knepley sp->functional[f].numQuadPoints = 1; 170529d8e7bcSMatthew G. Knepley ierr = PetscMalloc(sp->functional[f].numQuadPoints*dim * sizeof(PetscReal), &qpoints);CHKERRQ(ierr); 170629d8e7bcSMatthew G. Knepley ierr = PetscMalloc(sp->functional[f].numQuadPoints * sizeof(PetscReal), &qweights);CHKERRQ(ierr); 170759804f93SMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 170829d8e7bcSMatthew G. Knepley ierr = PetscSectionGetDof(csection, p, &dof);CHKERRQ(ierr); 170929d8e7bcSMatthew G. Knepley ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr); 171029d8e7bcSMatthew 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); 171129d8e7bcSMatthew G. Knepley for (d = 0; d < dof; ++d) {qpoints[d] = coords[off+d];} 171229d8e7bcSMatthew G. Knepley qweights[0] = 1.0; 171329d8e7bcSMatthew G. Knepley sp->functional[f].quadPoints = qpoints; 171429d8e7bcSMatthew G. Knepley sp->functional[f].quadWeights = qweights; 171529d8e7bcSMatthew G. Knepley ++f; 171659804f93SMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 171759804f93SMatthew G. Knepley } else if ((p >= pStart[1]) && (p < pEnd[1])) { 171859804f93SMatthew G. Knepley /* Edges */ 1719*d9dbb98eSMatthew G. Knepley PetscScalar *coords; 172059804f93SMatthew G. Knepley PetscInt k; 172159804f93SMatthew G. Knepley 172259804f93SMatthew G. Knepley if (order < 2) continue; 172359804f93SMatthew G. Knepley coords = NULL; 172459804f93SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, csection, coordinates, p, &n, &coords);CHKERRQ(ierr); 172559804f93SMatthew G. Knepley if (n != dim*2) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d has %d coordinate values instead of %d", p, n, dim*2); 172659804f93SMatthew G. Knepley for (k = 1; k < order; ++k) { 172759804f93SMatthew G. Knepley sp->functional[f].numQuadPoints = 1; 172859804f93SMatthew G. Knepley ierr = PetscMalloc(sp->functional[f].numQuadPoints*dim * sizeof(PetscReal), &qpoints);CHKERRQ(ierr); 172959804f93SMatthew G. Knepley ierr = PetscMalloc(sp->functional[f].numQuadPoints * sizeof(PetscReal), &qweights);CHKERRQ(ierr); 173059804f93SMatthew G. Knepley for (d = 0; d < dim; ++d) {qpoints[d] = k*(coords[1*dim+d] - coords[0*dim+d])/order + coords[0*dim+d];} 173159804f93SMatthew G. Knepley qweights[0] = 1.0; 173259804f93SMatthew G. Knepley sp->functional[f].quadPoints = qpoints; 173359804f93SMatthew G. Knepley sp->functional[f].quadWeights = qweights; 173459804f93SMatthew G. Knepley ++f; 173559804f93SMatthew G. Knepley } 173659804f93SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, csection, coordinates, p, &n, &coords);CHKERRQ(ierr); 173759804f93SMatthew G. Knepley } else if ((p >= pStart[depth-1]) && (p < pEnd[depth-1])) { 173859804f93SMatthew G. Knepley /* Faces */ 173959804f93SMatthew G. Knepley 174059804f93SMatthew G. Knepley if (order < 3) continue; 174159804f93SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to implement faces"); 174259804f93SMatthew G. Knepley } else if ((p >= pStart[depth]) && (p < pEnd[depth])) { 174359804f93SMatthew G. Knepley /* Cells */ 174459804f93SMatthew G. Knepley 174559804f93SMatthew G. Knepley if ((order > 0) && (order < 3)) continue; 174659804f93SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to implement cells"); 174729d8e7bcSMatthew G. Knepley } 174829d8e7bcSMatthew G. Knepley } 174959804f93SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, pStart[depth], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 175029d8e7bcSMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle cells with cone size %d", coneSize); 175159804f93SMatthew G. Knepley ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 175259804f93SMatthew G. Knepley if (f != pdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vector %d not equal to dimension %d", f, pdim); 175329d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 175429d8e7bcSMatthew G. Knepley } 175529d8e7bcSMatthew G. Knepley 175629d8e7bcSMatthew G. Knepley #undef __FUNCT__ 175729d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceGetDimension_Lagrange" 175829d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim) 175929d8e7bcSMatthew G. Knepley { 176029d8e7bcSMatthew G. Knepley PetscInt deg = sp->order; 176129d8e7bcSMatthew G. Knepley PetscReal D = 1.0; 176229d8e7bcSMatthew G. Knepley PetscInt n, i; 176329d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 176429d8e7bcSMatthew G. Knepley 176529d8e7bcSMatthew G. Knepley PetscFunctionBegin; 176629d8e7bcSMatthew G. Knepley /* TODO: Assumes simplices */ 176729d8e7bcSMatthew G. Knepley ierr = DMPlexGetDimension(sp->dm, &n);CHKERRQ(ierr); 176829d8e7bcSMatthew G. Knepley for (i = 1; i <= n; ++i) { 176929d8e7bcSMatthew G. Knepley D *= ((PetscReal) (deg+i))/i; 177029d8e7bcSMatthew G. Knepley } 177129d8e7bcSMatthew G. Knepley *dim = (PetscInt) (D + 0.5); 177229d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 177329d8e7bcSMatthew G. Knepley } 177429d8e7bcSMatthew G. Knepley 177529d8e7bcSMatthew G. Knepley #undef __FUNCT__ 177629d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceInitialize_Lagrange" 177729d8e7bcSMatthew G. Knepley PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 177829d8e7bcSMatthew G. Knepley { 177929d8e7bcSMatthew G. Knepley PetscFunctionBegin; 178029d8e7bcSMatthew G. Knepley sp->ops->setfromoptions = NULL; 178129d8e7bcSMatthew G. Knepley sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 178229d8e7bcSMatthew G. Knepley sp->ops->view = NULL; 178329d8e7bcSMatthew G. Knepley sp->ops->destroy = NULL; 178429d8e7bcSMatthew G. Knepley sp->ops->getdimension = PetscDualSpaceGetDimension_Lagrange; 178529d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 178629d8e7bcSMatthew G. Knepley } 178729d8e7bcSMatthew G. Knepley 178829d8e7bcSMatthew G. Knepley /*MC 178929d8e7bcSMatthew G. Knepley PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals 179029d8e7bcSMatthew G. Knepley 179129d8e7bcSMatthew G. Knepley Level: intermediate 179229d8e7bcSMatthew G. Knepley 179329d8e7bcSMatthew G. Knepley .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 179429d8e7bcSMatthew G. Knepley M*/ 179529d8e7bcSMatthew G. Knepley 179629d8e7bcSMatthew G. Knepley #undef __FUNCT__ 179729d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscDualSpaceCreate_Lagrange" 179829d8e7bcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 179929d8e7bcSMatthew G. Knepley { 180029d8e7bcSMatthew G. Knepley PetscDualSpace_Lag *lag; 180129d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 180229d8e7bcSMatthew G. Knepley 180329d8e7bcSMatthew G. Knepley PetscFunctionBegin; 180429d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 180529d8e7bcSMatthew G. Knepley ierr = PetscNewLog(sp, PetscDualSpace_Lag, &lag);CHKERRQ(ierr); 180629d8e7bcSMatthew G. Knepley sp->data = lag; 180729d8e7bcSMatthew G. Knepley 180829d8e7bcSMatthew G. Knepley /* lag->n = 0; */ 180929d8e7bcSMatthew G. Knepley 181029d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr); 181129d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 181229d8e7bcSMatthew G. Knepley } 181329d8e7bcSMatthew G. Knepley 181429d8e7bcSMatthew G. Knepley 181529d8e7bcSMatthew G. Knepley PetscInt PETSCFE_CLASSID = 0; 181629d8e7bcSMatthew G. Knepley 181729d8e7bcSMatthew G. Knepley #undef __FUNCT__ 181829d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEView" 181929d8e7bcSMatthew G. Knepley /*@C 182029d8e7bcSMatthew G. Knepley PetscFEView - Views a PetscFE 182129d8e7bcSMatthew G. Knepley 182229d8e7bcSMatthew G. Knepley Collective on PetscFE 182329d8e7bcSMatthew G. Knepley 182429d8e7bcSMatthew G. Knepley Input Parameter: 182529d8e7bcSMatthew G. Knepley + sp - the PetscFE object to view 182629d8e7bcSMatthew G. Knepley - v - the viewer 182729d8e7bcSMatthew G. Knepley 182829d8e7bcSMatthew G. Knepley Level: developer 182929d8e7bcSMatthew G. Knepley 183029d8e7bcSMatthew G. Knepley .seealso PetscFEDestroy() 183129d8e7bcSMatthew G. Knepley @*/ 183229d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEView(PetscFE fem, PetscViewer v) 183329d8e7bcSMatthew G. Knepley { 183429d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 183529d8e7bcSMatthew G. Knepley 183629d8e7bcSMatthew G. Knepley PetscFunctionBegin; 183729d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 183829d8e7bcSMatthew G. Knepley if (!v) { 183929d8e7bcSMatthew G. Knepley ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &v);CHKERRQ(ierr); 184029d8e7bcSMatthew G. Knepley } 184129d8e7bcSMatthew G. Knepley if (fem->ops->view) { 184229d8e7bcSMatthew G. Knepley ierr = (*fem->ops->view)(fem, v);CHKERRQ(ierr); 184329d8e7bcSMatthew G. Knepley } 184429d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 184529d8e7bcSMatthew G. Knepley } 184629d8e7bcSMatthew G. Knepley 184729d8e7bcSMatthew G. Knepley #undef __FUNCT__ 184829d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEDestroy" 184929d8e7bcSMatthew G. Knepley /*@ 185029d8e7bcSMatthew G. Knepley PetscFEDestroy - Destroys a PetscFE object 185129d8e7bcSMatthew G. Knepley 185229d8e7bcSMatthew G. Knepley Collective on PetscFE 185329d8e7bcSMatthew G. Knepley 185429d8e7bcSMatthew G. Knepley Input Parameter: 185529d8e7bcSMatthew G. Knepley . fem - the PetscFE object to destroy 185629d8e7bcSMatthew G. Knepley 185729d8e7bcSMatthew G. Knepley Level: developer 185829d8e7bcSMatthew G. Knepley 185929d8e7bcSMatthew G. Knepley .seealso PetscFEView() 186029d8e7bcSMatthew G. Knepley @*/ 186129d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEDestroy(PetscFE *fem) 186229d8e7bcSMatthew G. Knepley { 186329d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 186429d8e7bcSMatthew G. Knepley 186529d8e7bcSMatthew G. Knepley PetscFunctionBegin; 186629d8e7bcSMatthew G. Knepley if (!*fem) PetscFunctionReturn(0); 186729d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1); 186829d8e7bcSMatthew G. Knepley 186929d8e7bcSMatthew G. Knepley if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);} 187029d8e7bcSMatthew G. Knepley ((PetscObject) (*fem))->refct = 0; 187129d8e7bcSMatthew G. Knepley /* if memory was published with AMS then destroy it */ 187229d8e7bcSMatthew G. Knepley ierr = PetscObjectAMSViewOff((PetscObject) *fem);CHKERRQ(ierr); 187329d8e7bcSMatthew G. Knepley 187429d8e7bcSMatthew G. Knepley ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr); 187529d8e7bcSMatthew G. Knepley ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr); 187629d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 187729d8e7bcSMatthew G. Knepley } 187829d8e7bcSMatthew G. Knepley 187929d8e7bcSMatthew G. Knepley #undef __FUNCT__ 188029d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFECreate" 188129d8e7bcSMatthew G. Knepley /*@ 188229d8e7bcSMatthew G. Knepley PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType(). 188329d8e7bcSMatthew G. Knepley 188429d8e7bcSMatthew G. Knepley Collective on MPI_Comm 188529d8e7bcSMatthew G. Knepley 188629d8e7bcSMatthew G. Knepley Input Parameter: 188729d8e7bcSMatthew G. Knepley . comm - The communicator for the PetscFE object 188829d8e7bcSMatthew G. Knepley 188929d8e7bcSMatthew G. Knepley Output Parameter: 189029d8e7bcSMatthew G. Knepley . fem - The PetscFE object 189129d8e7bcSMatthew G. Knepley 189229d8e7bcSMatthew G. Knepley Level: beginner 189329d8e7bcSMatthew G. Knepley 189429d8e7bcSMatthew G. Knepley .seealso: PetscFESetType(), PETSCFEGALERKIN 189529d8e7bcSMatthew G. Knepley @*/ 189629d8e7bcSMatthew G. Knepley PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem) 189729d8e7bcSMatthew G. Knepley { 189829d8e7bcSMatthew G. Knepley PetscFE f; 189929d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 190029d8e7bcSMatthew G. Knepley 190129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 190229d8e7bcSMatthew G. Knepley PetscValidPointer(fem, 2); 190329d8e7bcSMatthew G. Knepley *fem = NULL; 190429d8e7bcSMatthew G. Knepley #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 190529d8e7bcSMatthew G. Knepley ierr = PetscFEInitializePackage();CHKERRQ(ierr); 190629d8e7bcSMatthew G. Knepley #endif 190729d8e7bcSMatthew G. Knepley 190829d8e7bcSMatthew G. Knepley ierr = PetscHeaderCreate(f, _p_PetscFE, struct _PetscFEOps, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr); 190929d8e7bcSMatthew G. Knepley ierr = PetscMemzero(f->ops, sizeof(struct _PetscFEOps));CHKERRQ(ierr); 191029d8e7bcSMatthew G. Knepley 191129d8e7bcSMatthew G. Knepley f->basisSpace = NULL; 191229d8e7bcSMatthew G. Knepley f->dualSpace = NULL; 1913*d9dbb98eSMatthew G. Knepley f->numComponents = 1; 191429d8e7bcSMatthew G. Knepley 191529d8e7bcSMatthew G. Knepley *fem = f; 191629d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 191729d8e7bcSMatthew G. Knepley } 191829d8e7bcSMatthew G. Knepley 191929d8e7bcSMatthew G. Knepley #undef __FUNCT__ 192029d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetDimension" 192129d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 192229d8e7bcSMatthew G. Knepley { 192329d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 192429d8e7bcSMatthew G. Knepley 192529d8e7bcSMatthew G. Knepley PetscFunctionBegin; 192629d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 192729d8e7bcSMatthew G. Knepley PetscValidPointer(dim, 2); 192829d8e7bcSMatthew G. Knepley ierr = PetscSpaceGetDimension(fem->basisSpace, dim);CHKERRQ(ierr); 192929d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 193029d8e7bcSMatthew G. Knepley } 193129d8e7bcSMatthew G. Knepley 193229d8e7bcSMatthew G. Knepley #undef __FUNCT__ 1933*d9dbb98eSMatthew G. Knepley #define __FUNCT__ "PetscFESetNumComponents" 1934*d9dbb98eSMatthew G. Knepley PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp) 1935*d9dbb98eSMatthew G. Knepley { 1936*d9dbb98eSMatthew G. Knepley PetscFunctionBegin; 1937*d9dbb98eSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1938*d9dbb98eSMatthew G. Knepley fem->numComponents = comp; 1939*d9dbb98eSMatthew G. Knepley PetscFunctionReturn(0); 1940*d9dbb98eSMatthew G. Knepley } 1941*d9dbb98eSMatthew G. Knepley 1942*d9dbb98eSMatthew G. Knepley #undef __FUNCT__ 194329d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetNumComponents" 194429d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp) 194529d8e7bcSMatthew G. Knepley { 194629d8e7bcSMatthew G. Knepley PetscFunctionBegin; 194729d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 194829d8e7bcSMatthew G. Knepley PetscValidPointer(comp, 2); 1949*d9dbb98eSMatthew G. Knepley *comp = fem->numComponents; 195029d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 195129d8e7bcSMatthew G. Knepley } 195229d8e7bcSMatthew G. Knepley 195329d8e7bcSMatthew G. Knepley #undef __FUNCT__ 195429d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetBasisSpace" 195529d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp) 195629d8e7bcSMatthew G. Knepley { 195729d8e7bcSMatthew G. Knepley PetscFunctionBegin; 195829d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 195929d8e7bcSMatthew G. Knepley PetscValidPointer(sp, 2); 196029d8e7bcSMatthew G. Knepley *sp = fem->basisSpace; 196129d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 196229d8e7bcSMatthew G. Knepley } 196329d8e7bcSMatthew G. Knepley 196429d8e7bcSMatthew G. Knepley #undef __FUNCT__ 196529d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFESetBasisSpace" 196629d8e7bcSMatthew G. Knepley PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp) 196729d8e7bcSMatthew G. Knepley { 196829d8e7bcSMatthew G. Knepley PetscFunctionBegin; 196929d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 197029d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2); 197129d8e7bcSMatthew G. Knepley fem->basisSpace = sp; 197229d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 197329d8e7bcSMatthew G. Knepley } 197429d8e7bcSMatthew G. Knepley 197529d8e7bcSMatthew G. Knepley #undef __FUNCT__ 197629d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetDualSpace" 197729d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp) 197829d8e7bcSMatthew G. Knepley { 197929d8e7bcSMatthew G. Knepley PetscFunctionBegin; 198029d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 198129d8e7bcSMatthew G. Knepley PetscValidPointer(sp, 2); 198229d8e7bcSMatthew G. Knepley *sp = fem->dualSpace; 198329d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 198429d8e7bcSMatthew G. Knepley } 198529d8e7bcSMatthew G. Knepley 198629d8e7bcSMatthew G. Knepley #undef __FUNCT__ 198729d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFESetDualSpace" 198829d8e7bcSMatthew G. Knepley PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp) 198929d8e7bcSMatthew G. Knepley { 199029d8e7bcSMatthew G. Knepley PetscFunctionBegin; 199129d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 199229d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 199329d8e7bcSMatthew G. Knepley fem->dualSpace = sp; 199429d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 199529d8e7bcSMatthew G. Knepley } 199629d8e7bcSMatthew G. Knepley 199729d8e7bcSMatthew G. Knepley #undef __FUNCT__ 199829d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFEGetTabulation" 199929d8e7bcSMatthew G. Knepley PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 200029d8e7bcSMatthew G. Knepley { 200129d8e7bcSMatthew G. Knepley DM dm; 200229d8e7bcSMatthew G. Knepley PetscReal *tmpB, *invV; 200329d8e7bcSMatthew G. Knepley PetscInt pdim; /* Dimension of FE space P */ 200429d8e7bcSMatthew G. Knepley PetscInt dim; /* Spatial dimension */ 2005*d9dbb98eSMatthew G. Knepley PetscInt comp; /* Field components */ 200629d8e7bcSMatthew G. Knepley PetscInt p, j, k; 200729d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 200829d8e7bcSMatthew G. Knepley 200929d8e7bcSMatthew G. Knepley PetscFunctionBegin; 201029d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 201129d8e7bcSMatthew G. Knepley PetscValidPointer(points, 3); 201229d8e7bcSMatthew G. Knepley if (B) PetscValidPointer(B, 4); 201329d8e7bcSMatthew G. Knepley if (D) PetscValidPointer(D, 5); 201429d8e7bcSMatthew G. Knepley if (H) PetscValidPointer(H, 6); 201529d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 201629d8e7bcSMatthew G. Knepley 201729d8e7bcSMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 201829d8e7bcSMatthew G. Knepley ierr = PetscSpaceGetDimension(fem->basisSpace, &pdim);CHKERRQ(ierr); 2019*d9dbb98eSMatthew G. Knepley ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr); 202029d8e7bcSMatthew 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); */ 202129d8e7bcSMatthew G. Knepley 202229d8e7bcSMatthew G. Knepley if (B) { 2023*d9dbb98eSMatthew G. Knepley ierr = DMGetWorkArray(dm, npoints*pdim*comp, PETSC_REAL, B);CHKERRQ(ierr); 202429d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr); 202529d8e7bcSMatthew G. Knepley } 202629d8e7bcSMatthew G. Knepley if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*dim, PETSC_REAL, D);CHKERRQ(ierr);} 202729d8e7bcSMatthew G. Knepley if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*dim*dim, PETSC_REAL, H);CHKERRQ(ierr);} 202829d8e7bcSMatthew G. Knepley ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr); 202929d8e7bcSMatthew G. Knepley 203029d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr); 203129d8e7bcSMatthew G. Knepley for (j = 0; j < pdim; ++j) { 203229d8e7bcSMatthew G. Knepley PetscReal *Bf; 203329d8e7bcSMatthew G. Knepley PetscQuadrature f; 203429d8e7bcSMatthew G. Knepley PetscInt q; 203529d8e7bcSMatthew G. Knepley 203629d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f); 203729d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr); 203829d8e7bcSMatthew G. Knepley ierr = PetscSpaceEvaluate(fem->basisSpace, f.numQuadPoints, f.quadPoints, Bf, NULL, NULL);CHKERRQ(ierr); 203929d8e7bcSMatthew G. Knepley for (k = 0; k < pdim; ++k) { 204029d8e7bcSMatthew G. Knepley /* n_j \cdot \phi_k */ 204129d8e7bcSMatthew G. Knepley invV[j*pdim+k] = 0.0; 204229d8e7bcSMatthew G. Knepley for (q = 0; q < f.numQuadPoints; ++q) { 204329d8e7bcSMatthew G. Knepley invV[j*pdim+k] += Bf[q*pdim+k]*f.quadWeights[q]; 204429d8e7bcSMatthew G. Knepley } 204529d8e7bcSMatthew G. Knepley } 204629d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr); 204729d8e7bcSMatthew G. Knepley } 204829d8e7bcSMatthew G. Knepley { 204929d8e7bcSMatthew G. Knepley PetscReal *work; 205029d8e7bcSMatthew G. Knepley PetscBLASInt *pivots; 205129d8e7bcSMatthew G. Knepley PetscBLASInt n = pdim, info; 205229d8e7bcSMatthew G. Knepley 205329d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr); 205429d8e7bcSMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr); 205529d8e7bcSMatthew G. Knepley PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invV, &n, pivots, &info)); 205629d8e7bcSMatthew G. Knepley PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invV, &n, pivots, work, &n, &info)); 205729d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr); 205829d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr); 205929d8e7bcSMatthew G. Knepley } 206029d8e7bcSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 206129d8e7bcSMatthew G. Knepley if (B) { 206229d8e7bcSMatthew G. Knepley /* Multiply by V^{-1} (pdim x pdim) */ 206329d8e7bcSMatthew G. Knepley for (j = 0; j < pdim; ++j) { 2064*d9dbb98eSMatthew G. Knepley const PetscInt i = (p*pdim + j)*comp; 2065*d9dbb98eSMatthew G. Knepley PetscInt c; 206629d8e7bcSMatthew G. Knepley 206729d8e7bcSMatthew G. Knepley (*B)[i] = 0.0; 206829d8e7bcSMatthew G. Knepley for (k = 0; k < pdim; ++k) { 206929d8e7bcSMatthew G. Knepley (*B)[i] += invV[k*pdim+j] * tmpB[p*pdim + k]; 207029d8e7bcSMatthew G. Knepley } 2071*d9dbb98eSMatthew G. Knepley for (c = 1; c < comp; ++c) { 2072*d9dbb98eSMatthew G. Knepley (*B)[i+c] = (*B)[i]; 2073*d9dbb98eSMatthew G. Knepley } 207429d8e7bcSMatthew G. Knepley } 207529d8e7bcSMatthew G. Knepley } 207629d8e7bcSMatthew G. Knepley } 207729d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr); 207829d8e7bcSMatthew G. Knepley if (B) { 207929d8e7bcSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr); 208029d8e7bcSMatthew G. Knepley } 208129d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 208229d8e7bcSMatthew G. Knepley } 208329d8e7bcSMatthew G. Knepley 208429d8e7bcSMatthew G. Knepley #undef __FUNCT__ 208529d8e7bcSMatthew G. Knepley #define __FUNCT__ "PetscFERestoreTabulation" 208629d8e7bcSMatthew G. Knepley PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **D2) 208729d8e7bcSMatthew G. Knepley { 208829d8e7bcSMatthew G. Knepley DM dm; 208929d8e7bcSMatthew G. Knepley PetscErrorCode ierr; 209029d8e7bcSMatthew G. Knepley 209129d8e7bcSMatthew G. Knepley PetscFunctionBegin; 209229d8e7bcSMatthew G. Knepley PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 209329d8e7bcSMatthew G. Knepley ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 209429d8e7bcSMatthew G. Knepley if (B) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, B);CHKERRQ(ierr);} 209529d8e7bcSMatthew G. Knepley if (D) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D);CHKERRQ(ierr);} 209629d8e7bcSMatthew G. Knepley if (D2) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D2);CHKERRQ(ierr);} 209729d8e7bcSMatthew G. Knepley PetscFunctionReturn(0); 209829d8e7bcSMatthew G. Knepley } 2099