xref: /petsc/src/ts/utils/dmplexts.c (revision f7148743b32e6acdc79bfeef32b74f4bf9614a56)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
337c2070cSMatthew G. Knepley #include <petsc/private/snesimpl.h>
4924a1b8fSMatthew G. Knepley #include <petscds.h>
56dbbd306SMatthew G. Knepley #include <petscfv.h>
66dbbd306SMatthew G. Knepley 
7254c1ad2SMatthew G. Knepley #undef __FUNCT__
86da023fcSToby Isaac #define __FUNCT__ "DMTSConvertPlex"
96da023fcSToby Isaac static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy)
106da023fcSToby Isaac {
116da023fcSToby Isaac   PetscBool      isPlex;
126da023fcSToby Isaac   PetscErrorCode ierr;
136da023fcSToby Isaac 
146da023fcSToby Isaac   PetscFunctionBegin;
156da023fcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
166da023fcSToby Isaac   if (isPlex) {
176da023fcSToby Isaac     *plex = dm;
186da023fcSToby Isaac     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
19*f7148743SMatthew G. Knepley   } else {
20*f7148743SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
21*f7148743SMatthew G. Knepley     if (!*plex) {
226da023fcSToby Isaac       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
23*f7148743SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
246da023fcSToby Isaac       if (copy) {
256da023fcSToby Isaac         PetscInt    i;
266da023fcSToby Isaac         PetscObject obj;
276da023fcSToby Isaac         const char *comps[3] = {"A","dmAux","dmCh"};
286da023fcSToby Isaac 
296da023fcSToby Isaac         ierr = DMCopyDMTS(dm, *plex);CHKERRQ(ierr);
306da023fcSToby Isaac         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
316da023fcSToby Isaac         for (i = 0; i < 3; i++) {
326da023fcSToby Isaac           ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr);
336da023fcSToby Isaac           ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr);
346da023fcSToby Isaac         }
356da023fcSToby Isaac       }
36*f7148743SMatthew G. Knepley     } else {
37*f7148743SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
38*f7148743SMatthew G. Knepley     }
396da023fcSToby Isaac   }
406da023fcSToby Isaac   PetscFunctionReturn(0);
416da023fcSToby Isaac }
426da023fcSToby Isaac 
436da023fcSToby Isaac 
446da023fcSToby Isaac #undef __FUNCT__
45b6aca0f9SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometryFVM"
46c510411aSMatthew G. Knepley /*@
47b6aca0f9SMatthew G. Knepley   DMPlexTSGetGeometryFVM - Return precomputed geometric data
48c510411aSMatthew G. Knepley 
49c510411aSMatthew G. Knepley   Input Parameter:
50c510411aSMatthew G. Knepley . dm - The DM
51c510411aSMatthew G. Knepley 
52c510411aSMatthew G. Knepley   Output Parameters:
53c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry
54c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry
55c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
56c510411aSMatthew G. Knepley 
57c510411aSMatthew G. Knepley   Level: developer
58c510411aSMatthew G. Knepley 
59c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal()
60c510411aSMatthew G. Knepley @*/
61b6aca0f9SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
62a0ac79e7SMatthew G. Knepley {
63a0ac79e7SMatthew G. Knepley   DMTS           dmts;
641faf85eaSMatthew G. Knepley   PetscObject    obj;
656da023fcSToby Isaac   DM             plex;
66a0ac79e7SMatthew G. Knepley   PetscErrorCode ierr;
67a0ac79e7SMatthew G. Knepley 
68a0ac79e7SMatthew G. Knepley   PetscFunctionBegin;
69924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
706da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
716da023fcSToby Isaac   ierr = DMGetDMTS(plex, &dmts);CHKERRQ(ierr);
72b6aca0f9SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
731faf85eaSMatthew G. Knepley   if (!obj) {
741faf85eaSMatthew G. Knepley     Vec cellgeom, facegeom;
751faf85eaSMatthew G. Knepley 
766da023fcSToby Isaac     ierr = DMPlexComputeGeometryFVM(plex, &cellgeom, &facegeom);CHKERRQ(ierr);
77b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
78b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
791faf85eaSMatthew G. Knepley     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
801faf85eaSMatthew G. Knepley     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
811faf85eaSMatthew G. Knepley   }
82b6aca0f9SMatthew G. Knepley   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
83b6aca0f9SMatthew G. Knepley   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
846da023fcSToby Isaac   if (minRadius) {ierr = DMPlexGetMinRadius(plex, minRadius);CHKERRQ(ierr);}
85b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
86924a1b8fSMatthew G. Knepley   PetscFunctionReturn(0);
87924a1b8fSMatthew G. Knepley }
88924a1b8fSMatthew G. Knepley 
89924a1b8fSMatthew G. Knepley #undef __FUNCT__
90c49ccbb3SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM"
91c49ccbb3SMatthew G. Knepley /*@C
92c49ccbb3SMatthew G. Knepley   DMPlexTSGetGradientDM - Return gradient data layout
93c49ccbb3SMatthew G. Knepley 
94c49ccbb3SMatthew G. Knepley   Input Parameters:
95c49ccbb3SMatthew G. Knepley + dm - The DM
96c49ccbb3SMatthew G. Knepley - fv - The PetscFV
97c49ccbb3SMatthew G. Knepley 
98c49ccbb3SMatthew G. Knepley   Output Parameter:
99c49ccbb3SMatthew G. Knepley . dmGrad - The layout for gradient values
100c49ccbb3SMatthew G. Knepley 
101c49ccbb3SMatthew G. Knepley   Level: developer
102c49ccbb3SMatthew G. Knepley 
103c49ccbb3SMatthew G. Knepley .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
104c49ccbb3SMatthew G. Knepley @*/
105c49ccbb3SMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
106c49ccbb3SMatthew G. Knepley {
107c49ccbb3SMatthew G. Knepley   DMTS           dmts;
108c49ccbb3SMatthew G. Knepley   PetscObject    obj;
109c49ccbb3SMatthew G. Knepley   PetscBool      computeGradients;
1106da023fcSToby Isaac   DM             plex;
111c49ccbb3SMatthew G. Knepley   PetscErrorCode ierr;
112c49ccbb3SMatthew G. Knepley 
113c49ccbb3SMatthew G. Knepley   PetscFunctionBegin;
114c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
115c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
116c49ccbb3SMatthew G. Knepley   PetscValidPointer(dmGrad,3);
117c49ccbb3SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
118c49ccbb3SMatthew G. Knepley   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
1196da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1206da023fcSToby Isaac   ierr = DMGetDMTS(plex, &dmts);CHKERRQ(ierr);
121c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
122c49ccbb3SMatthew G. Knepley   if (!obj) {
123c49ccbb3SMatthew G. Knepley     DM  dmGrad;
12408449791SMatthew G. Knepley     Vec faceGeometry, cellGeometry;
125c49ccbb3SMatthew G. Knepley 
1266da023fcSToby Isaac     ierr = DMPlexTSGetGeometryFVM(plex, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
1276da023fcSToby Isaac     ierr = DMPlexComputeGradientFVM(plex, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
128c49ccbb3SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
129924a1b8fSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
130c49ccbb3SMatthew G. Knepley   }
131c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
132b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
133254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
134254c1ad2SMatthew G. Knepley }
135254c1ad2SMatthew G. Knepley 
136254c1ad2SMatthew G. Knepley #undef __FUNCT__
137924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
13808449791SMatthew G. Knepley /*@
13908449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
14008449791SMatthew G. Knepley 
14108449791SMatthew G. Knepley   Input Parameters:
14208449791SMatthew G. Knepley + dm - The mesh
14308449791SMatthew G. Knepley . t - The time
14408449791SMatthew G. Knepley . locX  - Local solution
14508449791SMatthew G. Knepley - user - The user context
14608449791SMatthew G. Knepley 
14708449791SMatthew G. Knepley   Output Parameter:
1483b16df42SMatthew G. Knepley . F  - Global output vector
14908449791SMatthew G. Knepley 
15008449791SMatthew G. Knepley   Level: developer
15108449791SMatthew G. Knepley 
15208449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
15308449791SMatthew G. Knepley @*/
1543b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
155254c1ad2SMatthew G. Knepley {
1563b16df42SMatthew G. Knepley   Vec            locF;
157c4d4a4f8SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
1586da023fcSToby Isaac   DM             plex;
159254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
160254c1ad2SMatthew G. Knepley 
161254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
1626da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1636da023fcSToby Isaac   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
1646da023fcSToby Isaac   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
165c4d4a4f8SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1666da023fcSToby Isaac   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
1673b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
1686da023fcSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, NULL, locF, user);CHKERRQ(ierr);
1696da023fcSToby Isaac   ierr = DMLocalToGlobalBegin(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
1706da023fcSToby Isaac   ierr = DMLocalToGlobalEnd(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
1716da023fcSToby Isaac   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
172b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
173254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
174254c1ad2SMatthew G. Knepley }
175254c1ad2SMatthew G. Knepley 
176254c1ad2SMatthew G. Knepley #undef __FUNCT__
177ef68eab9SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeBoundary"
178c5d70e09SMatthew G. Knepley /*@
179ef68eab9SMatthew G. Knepley   DMPlexTSComputeBoundary - Insert the essential boundary values for the local input X and/or its time derivative X_t using pointwise functions specified by the user
180c5d70e09SMatthew G. Knepley 
181c5d70e09SMatthew G. Knepley   Input Parameters:
182c5d70e09SMatthew G. Knepley + dm - The mesh
183a40652d4SToby Isaac . t - The time
184a40652d4SToby Isaac . locX  - Local solution
185a40652d4SToby Isaac . locX_t - Local solution time derivative, or NULL
186c5d70e09SMatthew G. Knepley - user - The user context
187c5d70e09SMatthew G. Knepley 
188c5d70e09SMatthew G. Knepley   Level: developer
189c5d70e09SMatthew G. Knepley 
190a40652d4SToby Isaac .seealso: DMPlexComputeJacobianActionFEM()
191c5d70e09SMatthew G. Knepley @*/
192ef68eab9SMatthew G. Knepley PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
193c5d70e09SMatthew G. Knepley {
194c5d70e09SMatthew G. Knepley   DM             plex;
195ef68eab9SMatthew G. Knepley   Vec            faceGeometryFVM = NULL;
196ef68eab9SMatthew G. Knepley   PetscInt       Nf, f;
197c5d70e09SMatthew G. Knepley   PetscErrorCode ierr;
198c5d70e09SMatthew G. Knepley 
199c5d70e09SMatthew G. Knepley   PetscFunctionBegin;
200c5d70e09SMatthew G. Knepley   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
201ef68eab9SMatthew G. Knepley   ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr);
2029586001cSMatthew G. Knepley   if (!locX_t) {
2039586001cSMatthew G. Knepley     /* This is the RHS part */
204ef68eab9SMatthew G. Knepley     for (f = 0; f < Nf; f++) {
205ef68eab9SMatthew G. Knepley       PetscObject  obj;
206ef68eab9SMatthew G. Knepley       PetscClassId id;
207ef68eab9SMatthew G. Knepley 
208ef68eab9SMatthew G. Knepley       ierr = DMGetField(plex, f, &obj);CHKERRQ(ierr);
209ef68eab9SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
210ef68eab9SMatthew G. Knepley       if (id == PETSCFV_CLASSID) {
211ef68eab9SMatthew G. Knepley         ierr = DMPlexSNESGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
212ef68eab9SMatthew G. Knepley         break;
213ef68eab9SMatthew G. Knepley       }
214ef68eab9SMatthew G. Knepley     }
2159586001cSMatthew G. Knepley   }
2169586001cSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
217a40652d4SToby Isaac   /* TODO: locX_t */
218c5d70e09SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
219c5d70e09SMatthew G. Knepley   PetscFunctionReturn(0);
220c5d70e09SMatthew G. Knepley }
221c5d70e09SMatthew G. Knepley 
222c5d70e09SMatthew G. Knepley #undef __FUNCT__
22324cdb843SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
22424cdb843SMatthew G. Knepley /*@
22524cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
22624cdb843SMatthew G. Knepley 
22724cdb843SMatthew G. Knepley   Input Parameters:
22824cdb843SMatthew G. Knepley + dm - The mesh
22924cdb843SMatthew G. Knepley . t - The time
23008449791SMatthew G. Knepley . locX  - Local solution
23108449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
23224cdb843SMatthew G. Knepley - user - The user context
23324cdb843SMatthew G. Knepley 
23424cdb843SMatthew G. Knepley   Output Parameter:
23508449791SMatthew G. Knepley . locF  - Local output vector
23624cdb843SMatthew G. Knepley 
23724cdb843SMatthew G. Knepley   Level: developer
23824cdb843SMatthew G. Knepley 
23924cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
24024cdb843SMatthew G. Knepley @*/
24108449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
24224cdb843SMatthew G. Knepley {
243c4d4a4f8SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
2446da023fcSToby Isaac   DM             plex;
24524cdb843SMatthew G. Knepley   PetscErrorCode ierr;
24624cdb843SMatthew G. Knepley 
24724cdb843SMatthew G. Knepley   PetscFunctionBegin;
2486da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
2496da023fcSToby Isaac   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
2506da023fcSToby Isaac   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
251c4d4a4f8SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2526da023fcSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, locX_t, locF, user);CHKERRQ(ierr);
253b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
25424cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
25524cdb843SMatthew G. Knepley }
2567cdb2a12SMatthew G. Knepley 
2577cdb2a12SMatthew G. Knepley #undef __FUNCT__
258756a1f44SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIJacobianFEM"
259756a1f44SMatthew G. Knepley /*@
260756a1f44SMatthew G. Knepley   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
261756a1f44SMatthew G. Knepley 
262756a1f44SMatthew G. Knepley   Input Parameters:
263756a1f44SMatthew G. Knepley + dm - The mesh
264756a1f44SMatthew G. Knepley . t - The time
265756a1f44SMatthew G. Knepley . locX  - Local solution
266756a1f44SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
267756a1f44SMatthew G. Knepley . X_tshift - The multiplicative parameter for dF/du_t
268756a1f44SMatthew G. Knepley - user - The user context
269756a1f44SMatthew G. Knepley 
270756a1f44SMatthew G. Knepley   Output Parameter:
271756a1f44SMatthew G. Knepley . locF  - Local output vector
272756a1f44SMatthew G. Knepley 
273756a1f44SMatthew G. Knepley   Level: developer
274756a1f44SMatthew G. Knepley 
275756a1f44SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
276756a1f44SMatthew G. Knepley @*/
277756a1f44SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
278756a1f44SMatthew G. Knepley {
279756a1f44SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
280756a1f44SMatthew G. Knepley   DM             plex;
281756a1f44SMatthew G. Knepley   PetscErrorCode ierr;
282756a1f44SMatthew G. Knepley 
283756a1f44SMatthew G. Knepley   PetscFunctionBegin;
284756a1f44SMatthew G. Knepley   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
285756a1f44SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
286756a1f44SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
287756a1f44SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
288756a1f44SMatthew G. Knepley   ierr = DMPlexComputeJacobian_Internal(plex, cStart, cEnd, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr);
289756a1f44SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
290756a1f44SMatthew G. Knepley   PetscFunctionReturn(0);
291756a1f44SMatthew G. Knepley }
292756a1f44SMatthew G. Knepley 
293756a1f44SMatthew G. Knepley #undef __FUNCT__
2947cdb2a12SMatthew G. Knepley #define __FUNCT__ "DMTSCheckFromOptions"
2950163fd50SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
2967cdb2a12SMatthew G. Knepley {
2977cdb2a12SMatthew G. Knepley   DM             dm;
2987cdb2a12SMatthew G. Knepley   SNES           snes;
2991878804aSMatthew G. Knepley   Vec            sol;
3007cdb2a12SMatthew G. Knepley   PetscBool      check;
3017cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
3027cdb2a12SMatthew G. Knepley 
3037cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
304c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
3057cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
30672fd7fc8SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
30772fd7fc8SMatthew G. Knepley   ierr = TSSetSolution(ts, sol);CHKERRQ(ierr);
3087cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
3097cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
3107cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
3117cdb2a12SMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
3121878804aSMatthew G. Knepley   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
31372fd7fc8SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
3147cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
3157cdb2a12SMatthew G. Knepley }
316