xref: /petsc/src/ts/utils/dmplexts.c (revision a40652d484231dadb820fc3a9703f41e31810d43)
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);
196da023fcSToby Isaac   }
206da023fcSToby Isaac   else {
216da023fcSToby Isaac     ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
226da023fcSToby Isaac     if (copy) {
236da023fcSToby Isaac       PetscInt    i;
246da023fcSToby Isaac       PetscObject obj;
256da023fcSToby Isaac       const char *comps[3] = {"A","dmAux","dmCh"};
266da023fcSToby Isaac 
276da023fcSToby Isaac       ierr = DMCopyDMTS(dm,*plex);CHKERRQ(ierr);
286da023fcSToby Isaac       ierr = DMCopyDMSNES(dm,*plex);CHKERRQ(ierr);
296da023fcSToby Isaac       for (i = 0; i < 3; i++) {
306da023fcSToby Isaac         ierr = PetscObjectQuery((PetscObject)dm,comps[i],&obj);CHKERRQ(ierr);
316da023fcSToby Isaac         ierr = PetscObjectCompose((PetscObject)*plex,comps[i],obj);CHKERRQ(ierr);
326da023fcSToby Isaac       }
336da023fcSToby Isaac     }
346da023fcSToby Isaac   }
356da023fcSToby Isaac   PetscFunctionReturn(0);
366da023fcSToby Isaac }
376da023fcSToby Isaac 
386da023fcSToby Isaac 
396da023fcSToby Isaac #undef __FUNCT__
40b6aca0f9SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometryFVM"
41c510411aSMatthew G. Knepley /*@
42b6aca0f9SMatthew G. Knepley   DMPlexTSGetGeometryFVM - Return precomputed geometric data
43c510411aSMatthew G. Knepley 
44c510411aSMatthew G. Knepley   Input Parameter:
45c510411aSMatthew G. Knepley . dm - The DM
46c510411aSMatthew G. Knepley 
47c510411aSMatthew G. Knepley   Output Parameters:
48c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry
49c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry
50c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
51c510411aSMatthew G. Knepley 
52c510411aSMatthew G. Knepley   Level: developer
53c510411aSMatthew G. Knepley 
54c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal()
55c510411aSMatthew G. Knepley @*/
56b6aca0f9SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
57a0ac79e7SMatthew G. Knepley {
58a0ac79e7SMatthew G. Knepley   DMTS           dmts;
591faf85eaSMatthew G. Knepley   PetscObject    obj;
606da023fcSToby Isaac   DM             plex;
61a0ac79e7SMatthew G. Knepley   PetscErrorCode ierr;
62a0ac79e7SMatthew G. Knepley 
63a0ac79e7SMatthew G. Knepley   PetscFunctionBegin;
64924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
656da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
666da023fcSToby Isaac   ierr = DMGetDMTS(plex, &dmts);CHKERRQ(ierr);
67b6aca0f9SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
681faf85eaSMatthew G. Knepley   if (!obj) {
691faf85eaSMatthew G. Knepley     Vec cellgeom, facegeom;
701faf85eaSMatthew G. Knepley 
716da023fcSToby Isaac     ierr = DMPlexComputeGeometryFVM(plex, &cellgeom, &facegeom);CHKERRQ(ierr);
72b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
73b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
741faf85eaSMatthew G. Knepley     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
751faf85eaSMatthew G. Knepley     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
761faf85eaSMatthew G. Knepley   }
77b6aca0f9SMatthew G. Knepley   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
78b6aca0f9SMatthew G. Knepley   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
796da023fcSToby Isaac   if (minRadius) {ierr = DMPlexGetMinRadius(plex, minRadius);CHKERRQ(ierr);}
80b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
81924a1b8fSMatthew G. Knepley   PetscFunctionReturn(0);
82924a1b8fSMatthew G. Knepley }
83924a1b8fSMatthew G. Knepley 
84924a1b8fSMatthew G. Knepley #undef __FUNCT__
85c49ccbb3SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM"
86c49ccbb3SMatthew G. Knepley /*@C
87c49ccbb3SMatthew G. Knepley   DMPlexTSGetGradientDM - Return gradient data layout
88c49ccbb3SMatthew G. Knepley 
89c49ccbb3SMatthew G. Knepley   Input Parameters:
90c49ccbb3SMatthew G. Knepley + dm - The DM
91c49ccbb3SMatthew G. Knepley - fv - The PetscFV
92c49ccbb3SMatthew G. Knepley 
93c49ccbb3SMatthew G. Knepley   Output Parameter:
94c49ccbb3SMatthew G. Knepley . dmGrad - The layout for gradient values
95c49ccbb3SMatthew G. Knepley 
96c49ccbb3SMatthew G. Knepley   Level: developer
97c49ccbb3SMatthew G. Knepley 
98c49ccbb3SMatthew G. Knepley .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
99c49ccbb3SMatthew G. Knepley @*/
100c49ccbb3SMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
101c49ccbb3SMatthew G. Knepley {
102c49ccbb3SMatthew G. Knepley   DMTS           dmts;
103c49ccbb3SMatthew G. Knepley   PetscObject    obj;
104c49ccbb3SMatthew G. Knepley   PetscBool      computeGradients;
1056da023fcSToby Isaac   DM             plex;
106c49ccbb3SMatthew G. Knepley   PetscErrorCode ierr;
107c49ccbb3SMatthew G. Knepley 
108c49ccbb3SMatthew G. Knepley   PetscFunctionBegin;
109c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
110c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
111c49ccbb3SMatthew G. Knepley   PetscValidPointer(dmGrad,3);
112c49ccbb3SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
113c49ccbb3SMatthew G. Knepley   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
1146da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1156da023fcSToby Isaac   ierr = DMGetDMTS(plex, &dmts);CHKERRQ(ierr);
116c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
117c49ccbb3SMatthew G. Knepley   if (!obj) {
118c49ccbb3SMatthew G. Knepley     DM  dmGrad;
11908449791SMatthew G. Knepley     Vec faceGeometry, cellGeometry;
120c49ccbb3SMatthew G. Knepley 
1216da023fcSToby Isaac     ierr = DMPlexTSGetGeometryFVM(plex, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
1226da023fcSToby Isaac     ierr = DMPlexComputeGradientFVM(plex, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
123c49ccbb3SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
124924a1b8fSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
125c49ccbb3SMatthew G. Knepley   }
126c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
127b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
128254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
129254c1ad2SMatthew G. Knepley }
130254c1ad2SMatthew G. Knepley 
131254c1ad2SMatthew G. Knepley #undef __FUNCT__
132*a40652d4SToby Isaac #define __FUNCT__ "DMPlexTSComputeRHSBoundaryFVM"
133*a40652d4SToby Isaac /*@
134*a40652d4SToby Isaac   DMPlexTSComputeRHSBoundaryFVM - Insert the essential boundary conditions in the local input X using pointwise functions specified by the user
135*a40652d4SToby Isaac 
136*a40652d4SToby Isaac   Input Parameters:
137*a40652d4SToby Isaac + dm - The mesh
138*a40652d4SToby Isaac . t - The time
139*a40652d4SToby Isaac . locX  - Local solution
140*a40652d4SToby Isaac - user - The user context
141*a40652d4SToby Isaac 
142*a40652d4SToby Isaac   Level: developer
143*a40652d4SToby Isaac 
144*a40652d4SToby Isaac .seealso: DMPlexComputeJacobianActionFEM()
145*a40652d4SToby Isaac @*/
146*a40652d4SToby Isaac PetscErrorCode DMPlexTSComputeRHSBoundaryFVM(DM dm, PetscReal time, Vec locX, void *user)
147*a40652d4SToby Isaac {
148*a40652d4SToby Isaac   DM             plex;
149*a40652d4SToby Isaac   PetscErrorCode ierr;
150*a40652d4SToby Isaac 
151*a40652d4SToby Isaac   PetscFunctionBegin;
152*a40652d4SToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
153*a40652d4SToby Isaac   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, NULL, NULL, NULL);CHKERRQ(ierr);
154*a40652d4SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
155*a40652d4SToby Isaac   PetscFunctionReturn(0);
156*a40652d4SToby Isaac }
157*a40652d4SToby Isaac 
158*a40652d4SToby Isaac #undef __FUNCT__
159924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
16008449791SMatthew G. Knepley /*@
16108449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
16208449791SMatthew G. Knepley 
16308449791SMatthew G. Knepley   Input Parameters:
16408449791SMatthew G. Knepley + dm - The mesh
16508449791SMatthew G. Knepley . t - The time
16608449791SMatthew G. Knepley . locX  - Local solution
16708449791SMatthew G. Knepley - user - The user context
16808449791SMatthew G. Knepley 
16908449791SMatthew G. Knepley   Output Parameter:
1703b16df42SMatthew G. Knepley . F  - Global output vector
17108449791SMatthew G. Knepley 
17208449791SMatthew G. Knepley   Level: developer
17308449791SMatthew G. Knepley 
17408449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
17508449791SMatthew G. Knepley @*/
1763b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
177254c1ad2SMatthew G. Knepley {
1783b16df42SMatthew G. Knepley   Vec            locF;
179c4d4a4f8SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
1806da023fcSToby Isaac   DM             plex;
181254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
182254c1ad2SMatthew G. Knepley 
183254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
1846da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1856da023fcSToby Isaac   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
1866da023fcSToby Isaac   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
187c4d4a4f8SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1886da023fcSToby Isaac   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
1893b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
1906da023fcSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, NULL, locF, user);CHKERRQ(ierr);
1916da023fcSToby Isaac   ierr = DMLocalToGlobalBegin(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
1926da023fcSToby Isaac   ierr = DMLocalToGlobalEnd(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
1936da023fcSToby Isaac   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
194b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
195254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
196254c1ad2SMatthew G. Knepley }
197254c1ad2SMatthew G. Knepley 
198254c1ad2SMatthew G. Knepley #undef __FUNCT__
199*a40652d4SToby Isaac #define __FUNCT__ "DMPlexTSComputeIBoundaryFEM"
200*a40652d4SToby Isaac /*@
201*a40652d4SToby Isaac   DMPlexTSComputeIBoundaryFEM - Insert the essential boundary values for the local input X and/or its time derivative X_t using pointwise functions specified by the user
202*a40652d4SToby Isaac 
203*a40652d4SToby Isaac   Input Parameters:
204*a40652d4SToby Isaac + dm - The mesh
205*a40652d4SToby Isaac . t - The time
206*a40652d4SToby Isaac . locX  - Local solution
207*a40652d4SToby Isaac . locX_t - Local solution time derivative, or NULL
208*a40652d4SToby Isaac - user - The user context
209*a40652d4SToby Isaac 
210*a40652d4SToby Isaac   Level: developer
211*a40652d4SToby Isaac 
212*a40652d4SToby Isaac .seealso: DMPlexComputeJacobianActionFEM()
213*a40652d4SToby Isaac @*/
214*a40652d4SToby Isaac PetscErrorCode DMPlexTSComputeIBoundaryFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
215*a40652d4SToby Isaac {
216*a40652d4SToby Isaac   DM             plex;
217*a40652d4SToby Isaac   PetscErrorCode ierr;
218*a40652d4SToby Isaac 
219*a40652d4SToby Isaac   PetscFunctionBegin;
220*a40652d4SToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
221*a40652d4SToby Isaac   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, NULL, NULL, NULL);CHKERRQ(ierr);
222*a40652d4SToby Isaac   /* TODO: locX_t */
223*a40652d4SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
224*a40652d4SToby Isaac   PetscFunctionReturn(0);
225*a40652d4SToby Isaac }
226*a40652d4SToby Isaac 
227*a40652d4SToby Isaac #undef __FUNCT__
22824cdb843SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
22924cdb843SMatthew G. Knepley /*@
23024cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
23124cdb843SMatthew G. Knepley 
23224cdb843SMatthew G. Knepley   Input Parameters:
23324cdb843SMatthew G. Knepley + dm - The mesh
23424cdb843SMatthew G. Knepley . t - The time
23508449791SMatthew G. Knepley . locX  - Local solution
23608449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
23724cdb843SMatthew G. Knepley - user - The user context
23824cdb843SMatthew G. Knepley 
23924cdb843SMatthew G. Knepley   Output Parameter:
24008449791SMatthew G. Knepley . locF  - Local output vector
24124cdb843SMatthew G. Knepley 
24224cdb843SMatthew G. Knepley   Level: developer
24324cdb843SMatthew G. Knepley 
24424cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
24524cdb843SMatthew G. Knepley @*/
24608449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
24724cdb843SMatthew G. Knepley {
248c4d4a4f8SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
2496da023fcSToby Isaac   DM             plex;
25024cdb843SMatthew G. Knepley   PetscErrorCode ierr;
25124cdb843SMatthew G. Knepley 
25224cdb843SMatthew G. Knepley   PetscFunctionBegin;
2536da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
2546da023fcSToby Isaac   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
2556da023fcSToby Isaac   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
256c4d4a4f8SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2576da023fcSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, locX_t, locF, user);CHKERRQ(ierr);
258b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
25924cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
26024cdb843SMatthew G. Knepley }
2617cdb2a12SMatthew G. Knepley 
2627cdb2a12SMatthew G. Knepley #undef __FUNCT__
263756a1f44SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIJacobianFEM"
264756a1f44SMatthew G. Knepley /*@
265756a1f44SMatthew G. Knepley   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
266756a1f44SMatthew G. Knepley 
267756a1f44SMatthew G. Knepley   Input Parameters:
268756a1f44SMatthew G. Knepley + dm - The mesh
269756a1f44SMatthew G. Knepley . t - The time
270756a1f44SMatthew G. Knepley . locX  - Local solution
271756a1f44SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
272756a1f44SMatthew G. Knepley . X_tshift - The multiplicative parameter for dF/du_t
273756a1f44SMatthew G. Knepley - user - The user context
274756a1f44SMatthew G. Knepley 
275756a1f44SMatthew G. Knepley   Output Parameter:
276756a1f44SMatthew G. Knepley . locF  - Local output vector
277756a1f44SMatthew G. Knepley 
278756a1f44SMatthew G. Knepley   Level: developer
279756a1f44SMatthew G. Knepley 
280756a1f44SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
281756a1f44SMatthew G. Knepley @*/
282756a1f44SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
283756a1f44SMatthew G. Knepley {
284756a1f44SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
285756a1f44SMatthew G. Knepley   DM             plex;
286756a1f44SMatthew G. Knepley   PetscErrorCode ierr;
287756a1f44SMatthew G. Knepley 
288756a1f44SMatthew G. Knepley   PetscFunctionBegin;
289756a1f44SMatthew G. Knepley   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
290756a1f44SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
291756a1f44SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
292756a1f44SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
293756a1f44SMatthew G. Knepley   ierr = DMPlexComputeJacobian_Internal(plex, cStart, cEnd, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr);
294756a1f44SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
295756a1f44SMatthew G. Knepley   PetscFunctionReturn(0);
296756a1f44SMatthew G. Knepley }
297756a1f44SMatthew G. Knepley 
298756a1f44SMatthew G. Knepley #undef __FUNCT__
2997cdb2a12SMatthew G. Knepley #define __FUNCT__ "DMTSCheckFromOptions"
3000163fd50SMatthew 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)
3017cdb2a12SMatthew G. Knepley {
3027cdb2a12SMatthew G. Knepley   DM             dm;
3037cdb2a12SMatthew G. Knepley   SNES           snes;
3041878804aSMatthew G. Knepley   Vec            sol;
3057cdb2a12SMatthew G. Knepley   PetscBool      check;
3067cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
3077cdb2a12SMatthew G. Knepley 
3087cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
309c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
3107cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
31172fd7fc8SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
31272fd7fc8SMatthew G. Knepley   ierr = TSSetSolution(ts, sol);CHKERRQ(ierr);
3137cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
3147cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
3157cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
3167cdb2a12SMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
3171878804aSMatthew G. Knepley   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
31872fd7fc8SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
3197cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
3207cdb2a12SMatthew G. Knepley }
321