xref: /petsc/src/dm/dt/fe/impls/composite/fecomposite.c (revision 012bc3648b56b7befe35b618b2ec944650f76ae7)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/
320cf1dd8SToby Isaac #include <petscblaslapack.h>
4*012bc364SMatthew G. Knepley #include <petscdmplextransform.h>
520cf1dd8SToby Isaac 
62b99622eSMatthew G. Knepley static PetscErrorCode PetscFEDestroy_Composite(PetscFE fem)
720cf1dd8SToby Isaac {
820cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
920cf1dd8SToby Isaac   PetscErrorCode     ierr;
1020cf1dd8SToby Isaac 
1120cf1dd8SToby Isaac   PetscFunctionBegin;
1220cf1dd8SToby Isaac   ierr = PetscFree(cmp->embedding);CHKERRQ(ierr);
1320cf1dd8SToby Isaac   ierr = PetscFree(cmp);CHKERRQ(ierr);
1420cf1dd8SToby Isaac   PetscFunctionReturn(0);
1520cf1dd8SToby Isaac }
1620cf1dd8SToby Isaac 
172b99622eSMatthew G. Knepley static PetscErrorCode PetscFESetUp_Composite(PetscFE fem)
1820cf1dd8SToby Isaac {
1920cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
2020cf1dd8SToby Isaac   DM                 K;
21412e9a14SMatthew G. Knepley   DMPolytopeType     ct;
22*012bc364SMatthew G. Knepley   DMPlexTransform    tr;
2320cf1dd8SToby Isaac   PetscReal         *subpoint;
2420cf1dd8SToby Isaac   PetscBLASInt      *pivots;
2520cf1dd8SToby Isaac   PetscBLASInt       n, info;
2620cf1dd8SToby Isaac   PetscScalar       *work, *invVscalar;
2720cf1dd8SToby Isaac   PetscInt           dim, pdim, spdim, j, s;
28b4457527SToby Isaac   PetscSection       section;
2920cf1dd8SToby Isaac   PetscErrorCode     ierr;
3020cf1dd8SToby Isaac 
3120cf1dd8SToby Isaac   PetscFunctionBegin;
3220cf1dd8SToby Isaac   /* Get affine mapping from reference cell to each subcell */
3320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &K);CHKERRQ(ierr);
3420cf1dd8SToby Isaac   ierr = DMGetDimension(K, &dim);CHKERRQ(ierr);
35412e9a14SMatthew G. Knepley   ierr = DMPlexGetCellType(K, 0, &ct);CHKERRQ(ierr);
36*012bc364SMatthew G. Knepley   ierr = DMPlexTransformCreate(PETSC_COMM_SELF, &tr);CHKERRQ(ierr);
37*012bc364SMatthew G. Knepley   ierr = DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR);CHKERRQ(ierr);
38*012bc364SMatthew G. Knepley   ierr = DMPlexRefineRegularGetAffineTransforms(tr, ct, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);CHKERRQ(ierr);
39*012bc364SMatthew G. Knepley   ierr = DMPlexTransformDestroy(&tr);CHKERRQ(ierr);
4020cf1dd8SToby Isaac   /* Determine dof embedding into subelements */
4120cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
4220cf1dd8SToby Isaac   ierr = PetscSpaceGetDimension(fem->basisSpace, &spdim);CHKERRQ(ierr);
4320cf1dd8SToby Isaac   ierr = PetscMalloc1(cmp->numSubelements*spdim,&cmp->embedding);CHKERRQ(ierr);
4420cf1dd8SToby Isaac   ierr = DMGetWorkArray(K, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr);
45b4457527SToby Isaac   ierr = PetscDualSpaceGetSection(fem->dualSpace, &section);CHKERRQ(ierr);
4620cf1dd8SToby Isaac   for (s = 0; s < cmp->numSubelements; ++s) {
4720cf1dd8SToby Isaac     PetscInt sd = 0;
48b4457527SToby Isaac     PetscInt closureSize;
49b4457527SToby Isaac     PetscInt *closure = NULL;
5020cf1dd8SToby Isaac 
51b4457527SToby Isaac     ierr = DMPlexGetTransitiveClosure(K, s, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
52b4457527SToby Isaac     for (j = 0; j < closureSize; j++) {
53b4457527SToby Isaac       PetscInt point = closure[2*j];
54b4457527SToby Isaac       PetscInt dof, off, k;
5520cf1dd8SToby Isaac 
56b4457527SToby Isaac       ierr = PetscSectionGetDof(section, point, &dof);CHKERRQ(ierr);
57b4457527SToby Isaac       ierr = PetscSectionGetOffset(section, point, &off);CHKERRQ(ierr);
58b4457527SToby Isaac       for (k = 0; k < dof; k++) cmp->embedding[s*spdim+sd++] = off + k;
5920cf1dd8SToby Isaac     }
60b4457527SToby Isaac     ierr = DMPlexRestoreTransitiveClosure(K, s, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6120cf1dd8SToby Isaac     if (sd != spdim) SETERRQ3(PetscObjectComm((PetscObject) fem), PETSC_ERR_PLIB, "Subelement %d has %d dual basis vectors != %d", s, sd, spdim);
6220cf1dd8SToby Isaac   }
6320cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(K, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr);
6420cf1dd8SToby Isaac   /* Construct the change of basis from prime basis to nodal basis for each subelement */
6520cf1dd8SToby Isaac   ierr = PetscMalloc1(cmp->numSubelements*spdim*spdim,&fem->invV);CHKERRQ(ierr);
6620cf1dd8SToby Isaac   ierr = PetscMalloc2(spdim,&pivots,spdim,&work);CHKERRQ(ierr);
6720cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
6820cf1dd8SToby Isaac   ierr = PetscMalloc1(cmp->numSubelements*spdim*spdim,&invVscalar);CHKERRQ(ierr);
6920cf1dd8SToby Isaac #else
7020cf1dd8SToby Isaac   invVscalar = fem->invV;
7120cf1dd8SToby Isaac #endif
7220cf1dd8SToby Isaac   for (s = 0; s < cmp->numSubelements; ++s) {
7320cf1dd8SToby Isaac     for (j = 0; j < spdim; ++j) {
7420cf1dd8SToby Isaac       PetscReal       *Bf;
7520cf1dd8SToby Isaac       PetscQuadrature  f;
7620cf1dd8SToby Isaac       const PetscReal *points, *weights;
7720cf1dd8SToby Isaac       PetscInt         Nc, Nq, q, k;
7820cf1dd8SToby Isaac 
7920cf1dd8SToby Isaac       ierr = PetscDualSpaceGetFunctional(fem->dualSpace, cmp->embedding[s*spdim+j], &f);CHKERRQ(ierr);
8020cf1dd8SToby Isaac       ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
8120cf1dd8SToby Isaac       ierr = PetscMalloc1(f->numPoints*spdim*Nc,&Bf);CHKERRQ(ierr);
8220cf1dd8SToby Isaac       ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr);
8320cf1dd8SToby Isaac       for (k = 0; k < spdim; ++k) {
8420cf1dd8SToby Isaac         /* n_j \cdot \phi_k */
8520cf1dd8SToby Isaac         invVscalar[(s*spdim + j)*spdim+k] = 0.0;
8620cf1dd8SToby Isaac         for (q = 0; q < Nq; ++q) {
8720cf1dd8SToby Isaac           invVscalar[(s*spdim + j)*spdim+k] += Bf[q*spdim+k]*weights[q];
8820cf1dd8SToby Isaac         }
8920cf1dd8SToby Isaac       }
9020cf1dd8SToby Isaac       ierr = PetscFree(Bf);CHKERRQ(ierr);
9120cf1dd8SToby Isaac     }
9220cf1dd8SToby Isaac     n = spdim;
9320cf1dd8SToby Isaac     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, &invVscalar[s*spdim*spdim], &n, pivots, &info));
9420cf1dd8SToby Isaac     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, &invVscalar[s*spdim*spdim], &n, pivots, work, &n, &info));
9520cf1dd8SToby Isaac   }
9620cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
9720cf1dd8SToby Isaac   for (s = 0; s <cmp->numSubelements*spdim*spdim; s++) fem->invV[s] = PetscRealPart(invVscalar[s]);
9820cf1dd8SToby Isaac   ierr = PetscFree(invVscalar);CHKERRQ(ierr);
9920cf1dd8SToby Isaac #endif
10020cf1dd8SToby Isaac   ierr = PetscFree2(pivots,work);CHKERRQ(ierr);
10120cf1dd8SToby Isaac   PetscFunctionReturn(0);
10220cf1dd8SToby Isaac }
10320cf1dd8SToby Isaac 
104ef0bb6c7SMatthew G. Knepley static PetscErrorCode PetscFECreateTabulation_Composite(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
10520cf1dd8SToby Isaac {
10620cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
10720cf1dd8SToby Isaac   DM                 dm;
108412e9a14SMatthew G. Knepley   DMPolytopeType     ct;
10920cf1dd8SToby Isaac   PetscInt           pdim;  /* Dimension of FE space P */
11020cf1dd8SToby Isaac   PetscInt           spdim; /* Dimension of subelement FE space P */
11120cf1dd8SToby Isaac   PetscInt           dim;   /* Spatial dimension */
11220cf1dd8SToby Isaac   PetscInt           comp;  /* Field components */
11320cf1dd8SToby Isaac   PetscInt          *subpoints;
114ef0bb6c7SMatthew G. Knepley   PetscReal         *B = K >= 0 ? T->T[0] : NULL;
115ef0bb6c7SMatthew G. Knepley   PetscReal         *D = K >= 1 ? T->T[1] : NULL;
116ef0bb6c7SMatthew G. Knepley   PetscReal         *H = K >= 2 ? T->T[2] : NULL;
117ef0bb6c7SMatthew G. Knepley   PetscReal         *tmpB = NULL, *tmpD = NULL, *tmpH = NULL, *subpoint;
11820cf1dd8SToby Isaac   PetscInt           p, s, d, e, j, k;
11920cf1dd8SToby Isaac   PetscErrorCode     ierr;
12020cf1dd8SToby Isaac 
12120cf1dd8SToby Isaac   PetscFunctionBegin;
12220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
12320cf1dd8SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
124412e9a14SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, 0, &ct);CHKERRQ(ierr);
12520cf1dd8SToby Isaac   ierr = PetscSpaceGetDimension(fem->basisSpace, &spdim);CHKERRQ(ierr);
12620cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
12720cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr);
12820cf1dd8SToby Isaac   /* Divide points into subelements */
12920cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, npoints, MPIU_INT, &subpoints);CHKERRQ(ierr);
13020cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr);
13120cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
13220cf1dd8SToby Isaac     for (s = 0; s < cmp->numSubelements; ++s) {
13320cf1dd8SToby Isaac       PetscBool inside;
13420cf1dd8SToby Isaac 
13520cf1dd8SToby Isaac       /* Apply transform, and check that point is inside cell */
13620cf1dd8SToby Isaac       for (d = 0; d < dim; ++d) {
13720cf1dd8SToby Isaac         subpoint[d] = -1.0;
13820cf1dd8SToby 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]);
13920cf1dd8SToby Isaac       }
140*012bc364SMatthew G. Knepley       ierr = DMPolytopeInCellTest(ct, subpoint, &inside);CHKERRQ(ierr);
14120cf1dd8SToby Isaac       if (inside) {subpoints[p] = s; break;}
14220cf1dd8SToby Isaac     }
14320cf1dd8SToby Isaac     if (s >= cmp->numSubelements) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d was not found in any subelement", p);
14420cf1dd8SToby Isaac   }
14520cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr);
14620cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
147ef0bb6c7SMatthew G. Knepley   if (K >= 0) {ierr = DMGetWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
148ef0bb6c7SMatthew G. Knepley   if (K >= 1) {ierr = DMGetWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
149ef0bb6c7SMatthew G. Knepley   if (K >= 2) {ierr = DMGetWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
150ef0bb6c7SMatthew G. Knepley   ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH);CHKERRQ(ierr);
15120cf1dd8SToby Isaac   /* Translate to the nodal basis */
152ef0bb6c7SMatthew G. Knepley   if (K >= 0) {ierr = PetscArrayzero(B, npoints*pdim*comp);CHKERRQ(ierr);}
153ef0bb6c7SMatthew G. Knepley   if (K >= 1) {ierr = PetscArrayzero(D, npoints*pdim*comp*dim);CHKERRQ(ierr);}
154ef0bb6c7SMatthew G. Knepley   if (K >= 2) {ierr = PetscArrayzero(H, npoints*pdim*comp*dim*dim);CHKERRQ(ierr);}
15520cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
15620cf1dd8SToby Isaac     const PetscInt s = subpoints[p];
15720cf1dd8SToby Isaac 
15820cf1dd8SToby Isaac     if (B) {
15920cf1dd8SToby Isaac       /* Multiply by V^{-1} (spdim x spdim) */
16020cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
16120cf1dd8SToby Isaac         const PetscInt i = (p*pdim + cmp->embedding[s*spdim+j])*comp;
16220cf1dd8SToby Isaac 
16320cf1dd8SToby Isaac         B[i] = 0.0;
16420cf1dd8SToby Isaac         for (k = 0; k < spdim; ++k) {
16520cf1dd8SToby Isaac           B[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpB[p*spdim + k];
16620cf1dd8SToby Isaac         }
16720cf1dd8SToby Isaac       }
16820cf1dd8SToby Isaac     }
16920cf1dd8SToby Isaac     if (D) {
17020cf1dd8SToby Isaac       /* Multiply by V^{-1} (spdim x spdim) */
17120cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
17220cf1dd8SToby Isaac         for (d = 0; d < dim; ++d) {
17320cf1dd8SToby Isaac           const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim + d;
17420cf1dd8SToby Isaac 
17520cf1dd8SToby Isaac           D[i] = 0.0;
17620cf1dd8SToby Isaac           for (k = 0; k < spdim; ++k) {
17720cf1dd8SToby Isaac             D[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpD[(p*spdim + k)*dim + d];
17820cf1dd8SToby Isaac           }
17920cf1dd8SToby Isaac         }
18020cf1dd8SToby Isaac       }
18120cf1dd8SToby Isaac     }
18220cf1dd8SToby Isaac     if (H) {
18320cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
18420cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
18520cf1dd8SToby Isaac         for (d = 0; d < dim*dim; ++d) {
18620cf1dd8SToby Isaac           const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim*dim + d;
18720cf1dd8SToby Isaac 
18820cf1dd8SToby Isaac           H[i] = 0.0;
18920cf1dd8SToby Isaac           for (k = 0; k < spdim; ++k) {
19020cf1dd8SToby Isaac             H[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpH[(p*spdim + k)*dim*dim + d];
19120cf1dd8SToby Isaac           }
19220cf1dd8SToby Isaac         }
19320cf1dd8SToby Isaac       }
19420cf1dd8SToby Isaac     }
19520cf1dd8SToby Isaac   }
19620cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, npoints, MPIU_INT, &subpoints);CHKERRQ(ierr);
197ef0bb6c7SMatthew G. Knepley   if (K >= 0) {ierr = DMRestoreWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
198ef0bb6c7SMatthew G. Knepley   if (K >= 1) {ierr = DMRestoreWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
199ef0bb6c7SMatthew G. Knepley   if (K >= 2) {ierr = DMRestoreWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
20020cf1dd8SToby Isaac   PetscFunctionReturn(0);
20120cf1dd8SToby Isaac }
20220cf1dd8SToby Isaac 
2032b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Composite(PetscFE fem)
20420cf1dd8SToby Isaac {
20520cf1dd8SToby Isaac   PetscFunctionBegin;
20620cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
20720cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Composite;
20820cf1dd8SToby Isaac   fem->ops->view                    = NULL;
20920cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Composite;
21020cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
211ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Composite;
21220cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
21320cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
21420cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
21520cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
21620cf1dd8SToby Isaac   PetscFunctionReturn(0);
21720cf1dd8SToby Isaac }
21820cf1dd8SToby Isaac 
21920cf1dd8SToby Isaac /*MC
22020cf1dd8SToby Isaac   PETSCFECOMPOSITE = "composite" - A PetscFE object that represents a composite element
22120cf1dd8SToby Isaac 
22220cf1dd8SToby Isaac   Level: intermediate
22320cf1dd8SToby Isaac 
22420cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
22520cf1dd8SToby Isaac M*/
22620cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Composite(PetscFE fem)
22720cf1dd8SToby Isaac {
22820cf1dd8SToby Isaac   PetscFE_Composite *cmp;
22920cf1dd8SToby Isaac   PetscErrorCode     ierr;
23020cf1dd8SToby Isaac 
23120cf1dd8SToby Isaac   PetscFunctionBegin;
23220cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
23320cf1dd8SToby Isaac   ierr      = PetscNewLog(fem, &cmp);CHKERRQ(ierr);
23420cf1dd8SToby Isaac   fem->data = cmp;
23520cf1dd8SToby Isaac 
23620cf1dd8SToby Isaac   cmp->numSubelements = -1;
23720cf1dd8SToby Isaac   cmp->v0             = NULL;
23820cf1dd8SToby Isaac   cmp->jac            = NULL;
23920cf1dd8SToby Isaac 
24020cf1dd8SToby Isaac   ierr = PetscFEInitialize_Composite(fem);CHKERRQ(ierr);
24120cf1dd8SToby Isaac   PetscFunctionReturn(0);
24220cf1dd8SToby Isaac }
24320cf1dd8SToby Isaac 
24420cf1dd8SToby Isaac /*@C
24520cf1dd8SToby Isaac   PetscFECompositeGetMapping - Returns the mappings from the reference element to each subelement
24620cf1dd8SToby Isaac 
24720cf1dd8SToby Isaac   Not collective
24820cf1dd8SToby Isaac 
24920cf1dd8SToby Isaac   Input Parameter:
25020cf1dd8SToby Isaac . fem - The PetscFE object
25120cf1dd8SToby Isaac 
25220cf1dd8SToby Isaac   Output Parameters:
25320cf1dd8SToby Isaac + blockSize - The number of elements in a block
25420cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch
25520cf1dd8SToby Isaac . batchSize - The number of elements in a batch
25620cf1dd8SToby Isaac - numBatches - The number of batches in a chunk
25720cf1dd8SToby Isaac 
25820cf1dd8SToby Isaac   Level: intermediate
25920cf1dd8SToby Isaac 
26020cf1dd8SToby Isaac .seealso: PetscFECreate()
26120cf1dd8SToby Isaac @*/
26220cf1dd8SToby Isaac PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PetscInt *numSubelements, const PetscReal *v0[], const PetscReal *jac[], const PetscReal *invjac[])
26320cf1dd8SToby Isaac {
26420cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
26520cf1dd8SToby Isaac 
26620cf1dd8SToby Isaac   PetscFunctionBegin;
26720cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
26820cf1dd8SToby Isaac   if (numSubelements) {PetscValidPointer(numSubelements, 2); *numSubelements = cmp->numSubelements;}
26920cf1dd8SToby Isaac   if (v0)             {PetscValidPointer(v0, 3);             *v0             = cmp->v0;}
27020cf1dd8SToby Isaac   if (jac)            {PetscValidPointer(jac, 4);            *jac            = cmp->jac;}
27120cf1dd8SToby Isaac   if (invjac)         {PetscValidPointer(invjac, 5);         *invjac         = cmp->invjac;}
27220cf1dd8SToby Isaac   PetscFunctionReturn(0);
27320cf1dd8SToby Isaac }
274