xref: /petsc/src/dm/dt/fe/impls/composite/fecomposite.c (revision 6d8694c4fbab79f9439f1ad13c0386ba7ee1ca4b)
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>
4012bc364SMatthew G. Knepley #include <petscdmplextransform.h>
520cf1dd8SToby Isaac 
PetscFEDestroy_Composite(PetscFE fem)6d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEDestroy_Composite(PetscFE fem)
7d71ae5a4SJacob Faibussowitsch {
820cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *)fem->data;
920cf1dd8SToby Isaac 
1020cf1dd8SToby Isaac   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(PetscFree(cmp->embedding));
129566063dSJacob Faibussowitsch   PetscCall(PetscFree(cmp));
133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1420cf1dd8SToby Isaac }
1520cf1dd8SToby Isaac 
PetscFESetUp_Composite(PetscFE fem)16d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFESetUp_Composite(PetscFE fem)
17d71ae5a4SJacob Faibussowitsch {
1820cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *)fem->data;
1920cf1dd8SToby Isaac   DM                 K;
20412e9a14SMatthew G. Knepley   DMPolytopeType     ct;
21012bc364SMatthew G. Knepley   DMPlexTransform    tr;
2220cf1dd8SToby Isaac   PetscReal         *subpoint;
2320cf1dd8SToby Isaac   PetscBLASInt      *pivots;
2420cf1dd8SToby Isaac   PetscBLASInt       n, info;
2520cf1dd8SToby Isaac   PetscScalar       *work, *invVscalar;
2620cf1dd8SToby Isaac   PetscInt           dim, pdim, spdim, j, s;
27b4457527SToby Isaac   PetscSection       section;
2820cf1dd8SToby Isaac 
2920cf1dd8SToby Isaac   PetscFunctionBegin;
3020cf1dd8SToby Isaac   /* Get affine mapping from reference cell to each subcell */
319566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &K));
329566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(K, &dim));
339566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(K, 0, &ct));
349566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCreate(PETSC_COMM_SELF, &tr));
359566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR));
369566063dSJacob Faibussowitsch   PetscCall(DMPlexRefineRegularGetAffineTransforms(tr, ct, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac));
379566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformDestroy(&tr));
3820cf1dd8SToby Isaac   /* Determine dof embedding into subelements */
399566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
409566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(fem->basisSpace, &spdim));
419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cmp->numSubelements * spdim, &cmp->embedding));
429566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(K, dim, MPIU_REAL, &subpoint));
439566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetSection(fem->dualSpace, &section));
4420cf1dd8SToby Isaac   for (s = 0; s < cmp->numSubelements; ++s) {
4520cf1dd8SToby Isaac     PetscInt  sd = 0;
46b4457527SToby Isaac     PetscInt  closureSize;
47b4457527SToby Isaac     PetscInt *closure = NULL;
4820cf1dd8SToby Isaac 
499566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(K, s, PETSC_TRUE, &closureSize, &closure));
50b4457527SToby Isaac     for (j = 0; j < closureSize; j++) {
51b4457527SToby Isaac       PetscInt point = closure[2 * j];
52b4457527SToby Isaac       PetscInt dof, off, k;
5320cf1dd8SToby Isaac 
549566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, point, &dof));
559566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, point, &off));
56b4457527SToby Isaac       for (k = 0; k < dof; k++) cmp->embedding[s * spdim + sd++] = off + k;
5720cf1dd8SToby Isaac     }
589566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(K, s, PETSC_TRUE, &closureSize, &closure));
5963a3b9bcSJacob Faibussowitsch     PetscCheck(sd == spdim, PetscObjectComm((PetscObject)fem), PETSC_ERR_PLIB, "Subelement %" PetscInt_FMT " has %" PetscInt_FMT " dual basis vectors != %" PetscInt_FMT, s, sd, spdim);
6020cf1dd8SToby Isaac   }
619566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(K, dim, MPIU_REAL, &subpoint));
6220cf1dd8SToby Isaac   /* Construct the change of basis from prime basis to nodal basis for each subelement */
639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cmp->numSubelements * spdim * spdim, &fem->invV));
649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(spdim, &pivots, spdim, &work));
6520cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cmp->numSubelements * spdim * spdim, &invVscalar));
6720cf1dd8SToby Isaac #else
6820cf1dd8SToby Isaac   invVscalar = fem->invV;
6920cf1dd8SToby Isaac #endif
7020cf1dd8SToby Isaac   for (s = 0; s < cmp->numSubelements; ++s) {
7120cf1dd8SToby Isaac     for (j = 0; j < spdim; ++j) {
7220cf1dd8SToby Isaac       PetscReal       *Bf;
7320cf1dd8SToby Isaac       PetscQuadrature  f;
7420cf1dd8SToby Isaac       const PetscReal *points, *weights;
7520cf1dd8SToby Isaac       PetscInt         Nc, Nq, q, k;
7620cf1dd8SToby Isaac 
779566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, cmp->embedding[s * spdim + j], &f));
789566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights));
799566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(f->numPoints * spdim * Nc, &Bf));
809566063dSJacob Faibussowitsch       PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL));
8120cf1dd8SToby Isaac       for (k = 0; k < spdim; ++k) {
8220cf1dd8SToby Isaac         /* n_j \cdot \phi_k */
8320cf1dd8SToby Isaac         invVscalar[(s * spdim + j) * spdim + k] = 0.0;
84ad540459SPierre Jolivet         for (q = 0; q < Nq; ++q) invVscalar[(s * spdim + j) * spdim + k] += Bf[q * spdim + k] * weights[q];
8520cf1dd8SToby Isaac       }
869566063dSJacob Faibussowitsch       PetscCall(PetscFree(Bf));
8720cf1dd8SToby Isaac     }
886497c311SBarry Smith     PetscCall(PetscBLASIntCast(spdim, &n));
89792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, &invVscalar[s * spdim * spdim], &n, pivots, &info));
90792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&n, &invVscalar[s * spdim * spdim], &n, pivots, work, &n, &info));
9120cf1dd8SToby Isaac   }
9220cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
9320cf1dd8SToby Isaac   for (s = 0; s < cmp->numSubelements * spdim * spdim; s++) fem->invV[s] = PetscRealPart(invVscalar[s]);
949566063dSJacob Faibussowitsch   PetscCall(PetscFree(invVscalar));
9520cf1dd8SToby Isaac #endif
969566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pivots, work));
973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9820cf1dd8SToby Isaac }
9920cf1dd8SToby Isaac 
PetscFEComputeTabulation_Composite(PetscFE fem,PetscInt npoints,const PetscReal points[],PetscInt K,PetscTabulation T)100*ce78bad3SBarry Smith static PetscErrorCode PetscFEComputeTabulation_Composite(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
101d71ae5a4SJacob Faibussowitsch {
10220cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *)fem->data;
10320cf1dd8SToby Isaac   DM                 dm;
104412e9a14SMatthew G. Knepley   DMPolytopeType     ct;
10520cf1dd8SToby Isaac   PetscInt           pdim;  /* Dimension of FE space P */
10620cf1dd8SToby Isaac   PetscInt           spdim; /* Dimension of subelement FE space P */
10720cf1dd8SToby Isaac   PetscInt           dim;   /* Spatial dimension */
10820cf1dd8SToby Isaac   PetscInt           comp;  /* Field components */
10920cf1dd8SToby Isaac   PetscInt          *subpoints;
110ef0bb6c7SMatthew G. Knepley   PetscReal         *B    = K >= 0 ? T->T[0] : NULL;
111ef0bb6c7SMatthew G. Knepley   PetscReal         *D    = K >= 1 ? T->T[1] : NULL;
112ef0bb6c7SMatthew G. Knepley   PetscReal         *H    = K >= 2 ? T->T[2] : NULL;
113ef0bb6c7SMatthew G. Knepley   PetscReal         *tmpB = NULL, *tmpD = NULL, *tmpH = NULL, *subpoint;
11420cf1dd8SToby Isaac   PetscInt           p, s, d, e, j, k;
11520cf1dd8SToby Isaac 
11620cf1dd8SToby Isaac   PetscFunctionBegin;
1179566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
1189566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, 0, &ct));
1209566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(fem->basisSpace, &spdim));
1219566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
1229566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &comp));
12320cf1dd8SToby Isaac   /* Divide points into subelements */
1249566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, npoints, MPIU_INT, &subpoints));
1259566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, dim, MPIU_REAL, &subpoint));
12620cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
12720cf1dd8SToby Isaac     for (s = 0; s < cmp->numSubelements; ++s) {
12820cf1dd8SToby Isaac       PetscBool inside;
12920cf1dd8SToby Isaac 
13020cf1dd8SToby Isaac       /* Apply transform, and check that point is inside cell */
13120cf1dd8SToby Isaac       for (d = 0; d < dim; ++d) {
13220cf1dd8SToby Isaac         subpoint[d] = -1.0;
13320cf1dd8SToby 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]);
13420cf1dd8SToby Isaac       }
1359566063dSJacob Faibussowitsch       PetscCall(DMPolytopeInCellTest(ct, subpoint, &inside));
1369371c9d4SSatish Balay       if (inside) {
1379371c9d4SSatish Balay         subpoints[p] = s;
1389371c9d4SSatish Balay         break;
1399371c9d4SSatish Balay       }
14020cf1dd8SToby Isaac     }
14163a3b9bcSJacob Faibussowitsch     PetscCheck(s < cmp->numSubelements, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " was not found in any subelement", p);
14220cf1dd8SToby Isaac   }
1439566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, dim, MPIU_REAL, &subpoint));
14420cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
1459566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * spdim, MPIU_REAL, &tmpB));
1469566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * spdim * dim, MPIU_REAL, &tmpD));
1479566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * spdim * dim * dim, MPIU_REAL, &tmpH));
1489566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH));
14920cf1dd8SToby Isaac   /* Translate to the nodal basis */
1509566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(PetscArrayzero(B, npoints * pdim * comp));
1519566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(PetscArrayzero(D, npoints * pdim * comp * dim));
1529566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(PetscArrayzero(H, npoints * pdim * comp * dim * dim));
15320cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
15420cf1dd8SToby Isaac     const PetscInt s = subpoints[p];
15520cf1dd8SToby Isaac 
15620cf1dd8SToby Isaac     if (B) {
15720cf1dd8SToby Isaac       /* Multiply by V^{-1} (spdim x spdim) */
15820cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
15920cf1dd8SToby Isaac         const PetscInt i = (p * pdim + cmp->embedding[s * spdim + j]) * comp;
16020cf1dd8SToby Isaac 
16120cf1dd8SToby Isaac         B[i] = 0.0;
162ad540459SPierre Jolivet         for (k = 0; k < spdim; ++k) B[i] += fem->invV[(s * spdim + k) * spdim + j] * tmpB[p * spdim + k];
16320cf1dd8SToby Isaac       }
16420cf1dd8SToby Isaac     }
16520cf1dd8SToby Isaac     if (D) {
16620cf1dd8SToby Isaac       /* Multiply by V^{-1} (spdim x spdim) */
16720cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
16820cf1dd8SToby Isaac         for (d = 0; d < dim; ++d) {
16920cf1dd8SToby Isaac           const PetscInt i = ((p * pdim + cmp->embedding[s * spdim + j]) * comp + 0) * dim + d;
17020cf1dd8SToby Isaac 
17120cf1dd8SToby Isaac           D[i] = 0.0;
172ad540459SPierre Jolivet           for (k = 0; k < spdim; ++k) D[i] += fem->invV[(s * spdim + k) * spdim + j] * tmpD[(p * spdim + k) * dim + d];
17320cf1dd8SToby Isaac         }
17420cf1dd8SToby Isaac       }
17520cf1dd8SToby Isaac     }
17620cf1dd8SToby Isaac     if (H) {
17720cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
17820cf1dd8SToby Isaac       for (j = 0; j < spdim; ++j) {
17920cf1dd8SToby Isaac         for (d = 0; d < dim * dim; ++d) {
18020cf1dd8SToby Isaac           const PetscInt i = ((p * pdim + cmp->embedding[s * spdim + j]) * comp + 0) * dim * dim + d;
18120cf1dd8SToby Isaac 
18220cf1dd8SToby Isaac           H[i] = 0.0;
183ad540459SPierre Jolivet           for (k = 0; k < spdim; ++k) H[i] += fem->invV[(s * spdim + k) * spdim + j] * tmpH[(p * spdim + k) * dim * dim + d];
18420cf1dd8SToby Isaac         }
18520cf1dd8SToby Isaac       }
18620cf1dd8SToby Isaac     }
18720cf1dd8SToby Isaac   }
1889566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, npoints, MPIU_INT, &subpoints));
1899566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * spdim, MPIU_REAL, &tmpB));
1909566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * spdim * dim, MPIU_REAL, &tmpD));
1919566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * spdim * dim * dim, MPIU_REAL, &tmpH));
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19320cf1dd8SToby Isaac }
19420cf1dd8SToby Isaac 
PetscFEInitialize_Composite(PetscFE fem)195d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Composite(PetscFE fem)
196d71ae5a4SJacob Faibussowitsch {
19720cf1dd8SToby Isaac   PetscFunctionBegin;
19820cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
19920cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Composite;
20020cf1dd8SToby Isaac   fem->ops->view                    = NULL;
20120cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Composite;
20220cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
203*ce78bad3SBarry Smith   fem->ops->computetabulation       = PetscFEComputeTabulation_Composite;
20420cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
20520cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
20620cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */;
20720cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
2083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20920cf1dd8SToby Isaac }
21020cf1dd8SToby Isaac 
21120cf1dd8SToby Isaac /*MC
212dce8aebaSBarry Smith   PETSCFECOMPOSITE = "composite" - A `PetscFEType` that represents a composite element
21320cf1dd8SToby Isaac 
21420cf1dd8SToby Isaac   Level: intermediate
21520cf1dd8SToby Isaac 
216db781477SPatrick Sanan .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
21720cf1dd8SToby Isaac M*/
PetscFECreate_Composite(PetscFE fem)218d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Composite(PetscFE fem)
219d71ae5a4SJacob Faibussowitsch {
22020cf1dd8SToby Isaac   PetscFE_Composite *cmp;
22120cf1dd8SToby Isaac 
22220cf1dd8SToby Isaac   PetscFunctionBegin;
22320cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2244dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&cmp));
22520cf1dd8SToby Isaac   fem->data = cmp;
22620cf1dd8SToby Isaac 
22720cf1dd8SToby Isaac   cmp->numSubelements = -1;
22820cf1dd8SToby Isaac   cmp->v0             = NULL;
22920cf1dd8SToby Isaac   cmp->jac            = NULL;
23020cf1dd8SToby Isaac 
2319566063dSJacob Faibussowitsch   PetscCall(PetscFEInitialize_Composite(fem));
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23320cf1dd8SToby Isaac }
23420cf1dd8SToby Isaac 
23520cf1dd8SToby Isaac /*@C
23620cf1dd8SToby Isaac   PetscFECompositeGetMapping - Returns the mappings from the reference element to each subelement
23720cf1dd8SToby Isaac 
23820f4b53cSBarry Smith   Not Collective
23920cf1dd8SToby Isaac 
24020cf1dd8SToby Isaac   Input Parameter:
241dce8aebaSBarry Smith . fem - The `PetscFE` object
24220cf1dd8SToby Isaac 
24320cf1dd8SToby Isaac   Output Parameters:
244a4e35b19SJacob Faibussowitsch + numSubelements - The number of sub elements
245f13dfd9eSBarry Smith . v0             - The affine transformation for each element, an array of length $dim * Nc$. Pass `NULL` to ignore.
246f13dfd9eSBarry Smith . jac            - The Jacobian for each element, an array of length $dim^2 * Nc$. Pass `NULL` to ignore.
247f13dfd9eSBarry Smith - invjac         - The inverse of the Jacobian, an array of length $dim^2 * Nc$. Pass `NULL` to ignore.
24820cf1dd8SToby Isaac 
24920cf1dd8SToby Isaac   Level: intermediate
25020cf1dd8SToby Isaac 
251f13dfd9eSBarry Smith   Note:
252f13dfd9eSBarry Smith   Do not free the output arrays.
253f13dfd9eSBarry Smith 
254dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFECreate()`
25520cf1dd8SToby Isaac @*/
PetscFECompositeGetMapping(PetscFE fem,PeOp PetscInt * numSubelements,PeOp const PetscReal * v0[],PeOp const PetscReal * jac[],PeOp const PetscReal * invjac[])256*ce78bad3SBarry Smith PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PeOp PetscInt *numSubelements, PeOp const PetscReal *v0[], PeOp const PetscReal *jac[], PeOp const PetscReal *invjac[])
257d71ae5a4SJacob Faibussowitsch {
25820cf1dd8SToby Isaac   PetscFE_Composite *cmp = (PetscFE_Composite *)fem->data;
25920cf1dd8SToby Isaac 
26020cf1dd8SToby Isaac   PetscFunctionBegin;
26120cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2629371c9d4SSatish Balay   if (numSubelements) {
2634f572ea9SToby Isaac     PetscAssertPointer(numSubelements, 2);
2649371c9d4SSatish Balay     *numSubelements = cmp->numSubelements;
2659371c9d4SSatish Balay   }
2669371c9d4SSatish Balay   if (v0) {
2674f572ea9SToby Isaac     PetscAssertPointer(v0, 3);
2689371c9d4SSatish Balay     *v0 = cmp->v0;
2699371c9d4SSatish Balay   }
2709371c9d4SSatish Balay   if (jac) {
2714f572ea9SToby Isaac     PetscAssertPointer(jac, 4);
2729371c9d4SSatish Balay     *jac = cmp->jac;
2739371c9d4SSatish Balay   }
2749371c9d4SSatish Balay   if (invjac) {
2754f572ea9SToby Isaac     PetscAssertPointer(invjac, 5);
2769371c9d4SSatish Balay     *invjac = cmp->invjac;
2779371c9d4SSatish Balay   }
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27920cf1dd8SToby Isaac }
280