xref: /petsc/src/dm/dt/fe/impls/composite/fecomposite.c (revision 20cf1dd8cd656e27510885a71ea4be028bf48813)
1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2*20cf1dd8SToby Isaac #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/
3*20cf1dd8SToby Isaac #include <petsc/private/dmpleximpl.h> /* For CellRefiner */
4*20cf1dd8SToby Isaac #include <petscblaslapack.h>
5*20cf1dd8SToby Isaac 
6*20cf1dd8SToby Isaac PetscErrorCode PetscFEDestroy_Composite(PetscFE fem)
7*20cf1dd8SToby Isaac {
8*20cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
9*20cf1dd8SToby Isaac   PetscErrorCode     ierr;
10*20cf1dd8SToby Isaac 
11*20cf1dd8SToby Isaac   PetscFunctionBegin;
12*20cf1dd8SToby Isaac   ierr = CellRefinerRestoreAffineTransforms_Internal(cmp->cellRefiner, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);CHKERRQ(ierr);
13*20cf1dd8SToby Isaac   ierr = PetscFree(cmp->embedding);CHKERRQ(ierr);
14*20cf1dd8SToby Isaac   ierr = PetscFree(cmp);CHKERRQ(ierr);
15*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
16*20cf1dd8SToby Isaac }
17*20cf1dd8SToby Isaac 
18*20cf1dd8SToby Isaac PetscErrorCode PetscFESetUp_Composite(PetscFE fem)
19*20cf1dd8SToby Isaac {
20*20cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
21*20cf1dd8SToby Isaac   DM                 K;
22*20cf1dd8SToby Isaac   PetscReal         *subpoint;
23*20cf1dd8SToby Isaac   PetscBLASInt      *pivots;
24*20cf1dd8SToby Isaac   PetscBLASInt       n, info;
25*20cf1dd8SToby Isaac   PetscScalar       *work, *invVscalar;
26*20cf1dd8SToby Isaac   PetscInt           dim, pdim, spdim, j, s;
27*20cf1dd8SToby Isaac   PetscErrorCode     ierr;
28*20cf1dd8SToby Isaac 
29*20cf1dd8SToby Isaac   PetscFunctionBegin;
30*20cf1dd8SToby Isaac   /* Get affine mapping from reference cell to each subcell */
31*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &K);CHKERRQ(ierr);
32*20cf1dd8SToby Isaac   ierr = DMGetDimension(K, &dim);CHKERRQ(ierr);
33*20cf1dd8SToby Isaac   ierr = DMPlexGetCellRefiner_Internal(K, &cmp->cellRefiner);CHKERRQ(ierr);
34*20cf1dd8SToby Isaac   ierr = CellRefinerGetAffineTransforms_Internal(cmp->cellRefiner, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);CHKERRQ(ierr);
35*20cf1dd8SToby Isaac   /* Determine dof embedding into subelements */
36*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
37*20cf1dd8SToby Isaac   ierr = PetscSpaceGetDimension(fem->basisSpace, &spdim);CHKERRQ(ierr);
38*20cf1dd8SToby Isaac   ierr = PetscMalloc1(cmp->numSubelements*spdim,&cmp->embedding);CHKERRQ(ierr);
39*20cf1dd8SToby Isaac   ierr = DMGetWorkArray(K, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr);
40*20cf1dd8SToby Isaac   for (s = 0; s < cmp->numSubelements; ++s) {
41*20cf1dd8SToby Isaac     PetscInt sd = 0;
42*20cf1dd8SToby Isaac 
43*20cf1dd8SToby Isaac     for (j = 0; j < pdim; ++j) {
44*20cf1dd8SToby Isaac       PetscBool       inside;
45*20cf1dd8SToby Isaac       PetscQuadrature f;
46*20cf1dd8SToby Isaac       PetscInt        d, e;
47*20cf1dd8SToby Isaac 
48*20cf1dd8SToby Isaac       ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);CHKERRQ(ierr);
49*20cf1dd8SToby Isaac       /* Apply transform to first point, and check that point is inside subcell */
50*20cf1dd8SToby Isaac       for (d = 0; d < dim; ++d) {
51*20cf1dd8SToby Isaac         subpoint[d] = -1.0;
52*20cf1dd8SToby Isaac         for (e = 0; e < dim; ++e) subpoint[d] += cmp->invjac[(s*dim + d)*dim+e]*(f->points[e] - cmp->v0[s*dim+e]);
53*20cf1dd8SToby Isaac       }
54*20cf1dd8SToby Isaac       ierr = CellRefinerInCellTest_Internal(cmp->cellRefiner, subpoint, &inside);CHKERRQ(ierr);
55*20cf1dd8SToby Isaac       if (inside) {cmp->embedding[s*spdim+sd++] = j;}
56*20cf1dd8SToby Isaac     }
57*20cf1dd8SToby Isaac     if (sd != spdim) SETERRQ3(PetscObjectComm((PetscObject) fem), PETSC_ERR_PLIB, "Subelement %d has %d dual basis vectors != %d", s, sd, spdim);
58*20cf1dd8SToby Isaac   }
59*20cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(K, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr);
60*20cf1dd8SToby Isaac   /* Construct the change of basis from prime basis to nodal basis for each subelement */
61*20cf1dd8SToby Isaac   ierr = PetscMalloc1(cmp->numSubelements*spdim*spdim,&fem->invV);CHKERRQ(ierr);
62*20cf1dd8SToby Isaac   ierr = PetscMalloc2(spdim,&pivots,spdim,&work);CHKERRQ(ierr);
63*20cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
64*20cf1dd8SToby Isaac   ierr = PetscMalloc1(cmp->numSubelements*spdim*spdim,&invVscalar);CHKERRQ(ierr);
65*20cf1dd8SToby Isaac #else
66*20cf1dd8SToby Isaac   invVscalar = fem->invV;
67*20cf1dd8SToby Isaac #endif
68*20cf1dd8SToby Isaac   for (s = 0; s < cmp->numSubelements; ++s) {
69*20cf1dd8SToby Isaac     for (j = 0; j < spdim; ++j) {
70*20cf1dd8SToby Isaac       PetscReal       *Bf;
71*20cf1dd8SToby Isaac       PetscQuadrature  f;
72*20cf1dd8SToby Isaac       const PetscReal *points, *weights;
73*20cf1dd8SToby Isaac       PetscInt         Nc, Nq, q, k;
74*20cf1dd8SToby Isaac 
75*20cf1dd8SToby Isaac       ierr = PetscDualSpaceGetFunctional(fem->dualSpace, cmp->embedding[s*spdim+j], &f);CHKERRQ(ierr);
76*20cf1dd8SToby Isaac       ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
77*20cf1dd8SToby Isaac       ierr = PetscMalloc1(f->numPoints*spdim*Nc,&Bf);CHKERRQ(ierr);
78*20cf1dd8SToby Isaac       ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr);
79*20cf1dd8SToby Isaac       for (k = 0; k < spdim; ++k) {
80*20cf1dd8SToby Isaac         /* n_j \cdot \phi_k */
81*20cf1dd8SToby Isaac         invVscalar[(s*spdim + j)*spdim+k] = 0.0;
82*20cf1dd8SToby Isaac         for (q = 0; q < Nq; ++q) {
83*20cf1dd8SToby Isaac           invVscalar[(s*spdim + j)*spdim+k] += Bf[q*spdim+k]*weights[q];
84*20cf1dd8SToby Isaac         }
85*20cf1dd8SToby Isaac       }
86*20cf1dd8SToby Isaac       ierr = PetscFree(Bf);CHKERRQ(ierr);
87*20cf1dd8SToby Isaac     }
88*20cf1dd8SToby Isaac     n = spdim;
89*20cf1dd8SToby Isaac     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, &invVscalar[s*spdim*spdim], &n, pivots, &info));
90*20cf1dd8SToby Isaac     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, &invVscalar[s*spdim*spdim], &n, pivots, work, &n, &info));
91*20cf1dd8SToby Isaac   }
92*20cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
93*20cf1dd8SToby Isaac   for (s = 0; s <cmp->numSubelements*spdim*spdim; s++) fem->invV[s] = PetscRealPart(invVscalar[s]);
94*20cf1dd8SToby Isaac   ierr = PetscFree(invVscalar);CHKERRQ(ierr);
95*20cf1dd8SToby Isaac #endif
96*20cf1dd8SToby Isaac   ierr = PetscFree2(pivots,work);CHKERRQ(ierr);
97*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
98*20cf1dd8SToby Isaac }
99*20cf1dd8SToby Isaac 
100*20cf1dd8SToby Isaac PetscErrorCode PetscFEGetTabulation_Composite(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal *B, PetscReal *D, PetscReal *H)
101*20cf1dd8SToby Isaac {
102*20cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
103*20cf1dd8SToby Isaac   DM                 dm;
104*20cf1dd8SToby Isaac   PetscInt           pdim;  /* Dimension of FE space P */
105*20cf1dd8SToby Isaac   PetscInt           spdim; /* Dimension of subelement FE space P */
106*20cf1dd8SToby Isaac   PetscInt           dim;   /* Spatial dimension */
107*20cf1dd8SToby Isaac   PetscInt           comp;  /* Field components */
108*20cf1dd8SToby Isaac   PetscInt          *subpoints;
109*20cf1dd8SToby Isaac   PetscReal         *tmpB, *tmpD, *tmpH, *subpoint;
110*20cf1dd8SToby Isaac   PetscInt           p, s, d, e, j, k;
111*20cf1dd8SToby Isaac   PetscErrorCode     ierr;
112*20cf1dd8SToby Isaac 
113*20cf1dd8SToby Isaac   PetscFunctionBegin;
114*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
115*20cf1dd8SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
116*20cf1dd8SToby Isaac   ierr = PetscSpaceGetDimension(fem->basisSpace, &spdim);CHKERRQ(ierr);
117*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
118*20cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr);
119*20cf1dd8SToby Isaac   /* Divide points into subelements */
120*20cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, npoints, MPIU_INT, &subpoints);CHKERRQ(ierr);
121*20cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr);
122*20cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
123*20cf1dd8SToby Isaac     for (s = 0; s < cmp->numSubelements; ++s) {
124*20cf1dd8SToby Isaac       PetscBool inside;
125*20cf1dd8SToby Isaac 
126*20cf1dd8SToby Isaac       /* Apply transform, and check that point is inside cell */
127*20cf1dd8SToby Isaac       for (d = 0; d < dim; ++d) {
128*20cf1dd8SToby Isaac         subpoint[d] = -1.0;
129*20cf1dd8SToby Isaac         for (e = 0; e < dim; ++e) subpoint[d] += cmp->invjac[(s*dim + d)*dim+e]*(points[p*dim+e] - cmp->v0[s*dim+e]);
130*20cf1dd8SToby Isaac       }
131*20cf1dd8SToby Isaac       ierr = CellRefinerInCellTest_Internal(cmp->cellRefiner, subpoint, &inside);CHKERRQ(ierr);
132*20cf1dd8SToby Isaac       if (inside) {subpoints[p] = s; break;}
133*20cf1dd8SToby Isaac     }
134*20cf1dd8SToby Isaac     if (s >= cmp->numSubelements) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d was not found in any subelement", p);
135*20cf1dd8SToby Isaac   }
136*20cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr);
137*20cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
138*20cf1dd8SToby Isaac   if (B) {ierr = DMGetWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
139*20cf1dd8SToby Isaac   if (D) {ierr = DMGetWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
140*20cf1dd8SToby Isaac   if (H) {ierr = DMGetWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
141*20cf1dd8SToby Isaac   ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? tmpD : NULL, H ? tmpH : NULL);CHKERRQ(ierr);
142*20cf1dd8SToby Isaac   /* Translate to the nodal basis */
143*20cf1dd8SToby Isaac   if (B) {ierr = PetscMemzero(B, npoints*pdim*comp * sizeof(PetscReal));CHKERRQ(ierr);}
144*20cf1dd8SToby Isaac   if (D) {ierr = PetscMemzero(D, npoints*pdim*comp*dim * sizeof(PetscReal));CHKERRQ(ierr);}
145*20cf1dd8SToby Isaac   if (H) {ierr = PetscMemzero(H, npoints*pdim*comp*dim*dim * sizeof(PetscReal));CHKERRQ(ierr);}
146*20cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
147*20cf1dd8SToby Isaac     const PetscInt s = subpoints[p];
148*20cf1dd8SToby Isaac 
149*20cf1dd8SToby Isaac     if (B) {
150*20cf1dd8SToby Isaac       /* Multiply by V^{-1} (spdim x spdim) */
151*20cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
152*20cf1dd8SToby Isaac         const PetscInt i = (p*pdim + cmp->embedding[s*spdim+j])*comp;
153*20cf1dd8SToby Isaac 
154*20cf1dd8SToby Isaac         B[i] = 0.0;
155*20cf1dd8SToby Isaac         for (k = 0; k < spdim; ++k) {
156*20cf1dd8SToby Isaac           B[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpB[p*spdim + k];
157*20cf1dd8SToby Isaac         }
158*20cf1dd8SToby Isaac       }
159*20cf1dd8SToby Isaac     }
160*20cf1dd8SToby Isaac     if (D) {
161*20cf1dd8SToby Isaac       /* Multiply by V^{-1} (spdim x spdim) */
162*20cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
163*20cf1dd8SToby Isaac         for (d = 0; d < dim; ++d) {
164*20cf1dd8SToby Isaac           const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim + d;
165*20cf1dd8SToby Isaac 
166*20cf1dd8SToby Isaac           D[i] = 0.0;
167*20cf1dd8SToby Isaac           for (k = 0; k < spdim; ++k) {
168*20cf1dd8SToby Isaac             D[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpD[(p*spdim + k)*dim + d];
169*20cf1dd8SToby Isaac           }
170*20cf1dd8SToby Isaac         }
171*20cf1dd8SToby Isaac       }
172*20cf1dd8SToby Isaac     }
173*20cf1dd8SToby Isaac     if (H) {
174*20cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
175*20cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
176*20cf1dd8SToby Isaac         for (d = 0; d < dim*dim; ++d) {
177*20cf1dd8SToby Isaac           const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim*dim + d;
178*20cf1dd8SToby Isaac 
179*20cf1dd8SToby Isaac           H[i] = 0.0;
180*20cf1dd8SToby Isaac           for (k = 0; k < spdim; ++k) {
181*20cf1dd8SToby Isaac             H[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpH[(p*spdim + k)*dim*dim + d];
182*20cf1dd8SToby Isaac           }
183*20cf1dd8SToby Isaac         }
184*20cf1dd8SToby Isaac       }
185*20cf1dd8SToby Isaac     }
186*20cf1dd8SToby Isaac   }
187*20cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, npoints, MPIU_INT, &subpoints);CHKERRQ(ierr);
188*20cf1dd8SToby Isaac   if (B) {ierr = DMRestoreWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
189*20cf1dd8SToby Isaac   if (D) {ierr = DMRestoreWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
190*20cf1dd8SToby Isaac   if (H) {ierr = DMRestoreWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
191*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
192*20cf1dd8SToby Isaac }
193*20cf1dd8SToby Isaac 
194*20cf1dd8SToby Isaac PetscErrorCode PetscFEInitialize_Composite(PetscFE fem)
195*20cf1dd8SToby Isaac {
196*20cf1dd8SToby Isaac   PetscFunctionBegin;
197*20cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
198*20cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Composite;
199*20cf1dd8SToby Isaac   fem->ops->view                    = NULL;
200*20cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Composite;
201*20cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
202*20cf1dd8SToby Isaac   fem->ops->gettabulation           = PetscFEGetTabulation_Composite;
203*20cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
204*20cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
205*20cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
206*20cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
207*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
208*20cf1dd8SToby Isaac }
209*20cf1dd8SToby Isaac 
210*20cf1dd8SToby Isaac /*MC
211*20cf1dd8SToby Isaac   PETSCFECOMPOSITE = "composite" - A PetscFE object that represents a composite element
212*20cf1dd8SToby Isaac 
213*20cf1dd8SToby Isaac   Level: intermediate
214*20cf1dd8SToby Isaac 
215*20cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
216*20cf1dd8SToby Isaac M*/
217*20cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Composite(PetscFE fem)
218*20cf1dd8SToby Isaac {
219*20cf1dd8SToby Isaac   PetscFE_Composite *cmp;
220*20cf1dd8SToby Isaac   PetscErrorCode     ierr;
221*20cf1dd8SToby Isaac 
222*20cf1dd8SToby Isaac   PetscFunctionBegin;
223*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
224*20cf1dd8SToby Isaac   ierr      = PetscNewLog(fem, &cmp);CHKERRQ(ierr);
225*20cf1dd8SToby Isaac   fem->data = cmp;
226*20cf1dd8SToby Isaac 
227*20cf1dd8SToby Isaac   cmp->cellRefiner    = REFINER_NOOP;
228*20cf1dd8SToby Isaac   cmp->numSubelements = -1;
229*20cf1dd8SToby Isaac   cmp->v0             = NULL;
230*20cf1dd8SToby Isaac   cmp->jac            = NULL;
231*20cf1dd8SToby Isaac 
232*20cf1dd8SToby Isaac   ierr = PetscFEInitialize_Composite(fem);CHKERRQ(ierr);
233*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
234*20cf1dd8SToby Isaac }
235*20cf1dd8SToby Isaac 
236*20cf1dd8SToby Isaac /*@C
237*20cf1dd8SToby Isaac   PetscFECompositeGetMapping - Returns the mappings from the reference element to each subelement
238*20cf1dd8SToby Isaac 
239*20cf1dd8SToby Isaac   Not collective
240*20cf1dd8SToby Isaac 
241*20cf1dd8SToby Isaac   Input Parameter:
242*20cf1dd8SToby Isaac . fem - The PetscFE object
243*20cf1dd8SToby Isaac 
244*20cf1dd8SToby Isaac   Output Parameters:
245*20cf1dd8SToby Isaac + blockSize - The number of elements in a block
246*20cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch
247*20cf1dd8SToby Isaac . batchSize - The number of elements in a batch
248*20cf1dd8SToby Isaac - numBatches - The number of batches in a chunk
249*20cf1dd8SToby Isaac 
250*20cf1dd8SToby Isaac   Level: intermediate
251*20cf1dd8SToby Isaac 
252*20cf1dd8SToby Isaac .seealso: PetscFECreate()
253*20cf1dd8SToby Isaac @*/
254*20cf1dd8SToby Isaac PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PetscInt *numSubelements, const PetscReal *v0[], const PetscReal *jac[], const PetscReal *invjac[])
255*20cf1dd8SToby Isaac {
256*20cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
257*20cf1dd8SToby Isaac 
258*20cf1dd8SToby Isaac   PetscFunctionBegin;
259*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
260*20cf1dd8SToby Isaac   if (numSubelements) {PetscValidPointer(numSubelements, 2); *numSubelements = cmp->numSubelements;}
261*20cf1dd8SToby Isaac   if (v0)             {PetscValidPointer(v0, 3);             *v0             = cmp->v0;}
262*20cf1dd8SToby Isaac   if (jac)            {PetscValidPointer(jac, 4);            *jac            = cmp->jac;}
263*20cf1dd8SToby Isaac   if (invjac)         {PetscValidPointer(invjac, 5);         *invjac         = cmp->invjac;}
264*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
265*20cf1dd8SToby Isaac }
266