xref: /petsc/src/ts/utils/dmplexts.c (revision a2fddd78f770fa4fc19a8af67e65be331f27d92b)
1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2 #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
3 #include <petsc/private/snesimpl.h>
4 #include <petscds.h>
5 #include <petscfv.h>
6 
7 static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy)
8 {
9   PetscBool      isPlex;
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
14   if (isPlex) {
15     *plex = dm;
16     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
17   } else {
18     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
19     if (!*plex) {
20       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
21       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
22       if (copy) {
23         ierr = DMCopyDMTS(dm, *plex);CHKERRQ(ierr);
24         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
25         ierr = DMCopyAuxiliaryVec(dm, *plex);CHKERRQ(ierr);
26       }
27     } else {
28       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
29     }
30   }
31   PetscFunctionReturn(0);
32 }
33 
34 /*@
35   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
36 
37   Input Parameters:
38 + dm - The mesh
39 . t - The time
40 . locX  - Local solution
41 - user - The user context
42 
43   Output Parameter:
44 . F  - Global output vector
45 
46   Level: developer
47 
48 .seealso: DMPlexComputeJacobianActionFEM()
49 @*/
50 PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
51 {
52   Vec            locF;
53   IS             cellIS;
54   DM             plex;
55   PetscInt       depth;
56   PetscFormKey key = {NULL, 0, 0};
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
61   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
62   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
63   if (!cellIS) {
64     ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);
65   }
66   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
67   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
68   ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locF, user);CHKERRQ(ierr);
69   ierr = DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F);CHKERRQ(ierr);
70   ierr = DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F);CHKERRQ(ierr);
71   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
72   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
73   ierr = DMDestroy(&plex);CHKERRQ(ierr);
74   PetscFunctionReturn(0);
75 }
76 
77 /*@
78   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
79 
80   Input Parameters:
81 + dm - The mesh
82 . t - The time
83 . locX  - Local solution
84 . locX_t - Local solution time derivative, or NULL
85 - user - The user context
86 
87   Level: developer
88 
89 .seealso: DMPlexComputeJacobianActionFEM()
90 @*/
91 PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
92 {
93   DM             plex;
94   Vec            faceGeometryFVM = NULL;
95   PetscInt       Nf, f;
96   PetscErrorCode ierr;
97 
98   PetscFunctionBegin;
99   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
100   ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr);
101   if (!locX_t) {
102     /* This is the RHS part */
103     for (f = 0; f < Nf; f++) {
104       PetscObject  obj;
105       PetscClassId id;
106 
107       ierr = DMGetField(plex, f, NULL, &obj);CHKERRQ(ierr);
108       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
109       if (id == PETSCFV_CLASSID) {
110         ierr = DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
111         break;
112       }
113     }
114   }
115   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
116   ierr = DMPlexInsertTimeDerivativeBoundaryValues(plex, PETSC_TRUE, locX_t, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
117   ierr = DMDestroy(&plex);CHKERRQ(ierr);
118   PetscFunctionReturn(0);
119 }
120 
121 /*@
122   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
123 
124   Input Parameters:
125 + dm - The mesh
126 . t - The time
127 . locX  - Local solution
128 . locX_t - Local solution time derivative, or NULL
129 - user - The user context
130 
131   Output Parameter:
132 . locF  - Local output vector
133 
134   Level: developer
135 
136 .seealso: DMPlexComputeJacobianActionFEM()
137 @*/
138 PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
139 {
140   DM             plex;
141   IS             allcellIS;
142   PetscInt       Nds, s;
143   PetscErrorCode ierr;
144 
145   PetscFunctionBegin;
146   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
147   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
148   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
149   for (s = 0; s < Nds; ++s) {
150     PetscDS          ds;
151     IS               cellIS;
152     PetscFormKey key;
153 
154     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
155     key.value = 0;
156     key.field = 0;
157     key.part  = 0;
158     if (!key.label) {
159       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
160       cellIS = allcellIS;
161     } else {
162       IS pointIS;
163 
164       key.value = 1;
165       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
166       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
167       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
168     }
169     ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, locX_t, time, locF, user);CHKERRQ(ierr);
170     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
171   }
172   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
173   ierr = DMDestroy(&plex);CHKERRQ(ierr);
174   PetscFunctionReturn(0);
175 }
176 
177 /*@
178   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
179 
180   Input Parameters:
181 + dm - The mesh
182 . t - The time
183 . locX  - Local solution
184 . locX_t - Local solution time derivative, or NULL
185 . X_tshift - The multiplicative parameter for dF/du_t
186 - user - The user context
187 
188   Output Parameter:
189 . locF  - Local output vector
190 
191   Level: developer
192 
193 .seealso: DMPlexComputeJacobianActionFEM()
194 @*/
195 PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
196 {
197   DM             plex;
198   IS             allcellIS;
199   PetscBool      hasJac, hasPrec;
200   PetscInt       Nds, s;
201   PetscErrorCode ierr;
202 
203   PetscFunctionBegin;
204   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
205   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
206   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
207   for (s = 0; s < Nds; ++s) {
208     PetscDS          ds;
209     IS               cellIS;
210     PetscFormKey key;
211 
212     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
213     key.value = 0;
214     key.field = 0;
215     key.part  = 0;
216     if (!key.label) {
217       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
218       cellIS = allcellIS;
219     } else {
220       IS pointIS;
221 
222       key.value = 1;
223       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
224       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
225       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
226     }
227     if (!s) {
228       ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
229       ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
230       if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
231       ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
232     }
233     ierr = DMPlexComputeJacobian_Internal(plex, key, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr);
234     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
235   }
236   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
237   ierr = DMDestroy(&plex);CHKERRQ(ierr);
238   PetscFunctionReturn(0);
239 }
240 
241 /*@C
242   DMTSCheckResidual - Check the residual of the exact solution
243 
244   Input Parameters:
245 + ts  - the TS object
246 . dm  - the DM
247 . t   - the time
248 . u   - a DM vector
249 . u_t - a DM vector
250 - tol - A tolerance for the check, or -1 to print the results instead
251 
252   Output Parameters:
253 . residual - The residual norm of the exact solution, or NULL
254 
255   Level: developer
256 
257 .seealso: DNTSCheckFromOptions(), DMTSCheckJacobian(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
258 @*/
259 PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual)
260 {
261   MPI_Comm       comm;
262   Vec            r;
263   PetscReal      res;
264   PetscErrorCode ierr;
265 
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
268   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
269   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
270   if (residual) PetscValidRealPointer(residual, 7);
271   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
272   ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr);
273   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
274   ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr);
275   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
276   if (tol >= 0.0) {
277     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
278   } else if (residual) {
279     *residual = res;
280   } else {
281     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
282     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
283     ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", (PetscObject) dm);CHKERRQ(ierr);
284     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
285     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
286     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
287     ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", NULL);CHKERRQ(ierr);
288   }
289   ierr = VecDestroy(&r);CHKERRQ(ierr);
290   PetscFunctionReturn(0);
291 }
292 
293 /*@C
294   DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
295 
296   Input Parameters:
297 + ts  - the TS object
298 . dm  - the DM
299 . t   - the time
300 . u   - a DM vector
301 . u_t - a DM vector
302 - tol - A tolerance for the check, or -1 to print the results instead
303 
304   Output Parameters:
305 + isLinear - Flag indicaing that the function looks linear, or NULL
306 - convRate - The rate of convergence of the linear model, or NULL
307 
308   Level: developer
309 
310 .seealso: DNTSCheckFromOptions(), DMTSCheckResidual(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
311 @*/
312 PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
313 {
314   MPI_Comm       comm;
315   PetscDS        ds;
316   Mat            J, M;
317   MatNullSpace   nullspace;
318   PetscReal      dt, shift, slope, intercept;
319   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
320   PetscErrorCode ierr;
321 
322   PetscFunctionBegin;
323   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
324   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
325   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
326   if (isLinear) PetscValidBoolPointer(isLinear, 7);
327   if (convRate) PetscValidRealPointer(convRate, 8);
328   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
329   ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr);
330   /* Create and view matrices */
331   ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr);
332   shift = 1.0/dt;
333   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
334   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
335   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
336   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
337   if (hasJac && hasPrec) {
338     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
339     ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE);CHKERRQ(ierr);
340     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
341     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
342     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
343     ierr = MatDestroy(&M);CHKERRQ(ierr);
344   } else {
345     ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE);CHKERRQ(ierr);
346   }
347   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
348   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
349   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
350   /* Check nullspace */
351   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
352   if (nullspace) {
353     PetscBool isNull;
354     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
355     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
356   }
357   /* Taylor test */
358   {
359     PetscRandom rand;
360     Vec         du, uhat, uhat_t, r, rhat, df;
361     PetscReal   h;
362     PetscReal  *es, *hs, *errors;
363     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
364     PetscInt    Nv, v;
365 
366     /* Choose a perturbation direction */
367     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
368     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
369     ierr = VecSetRandom(du, rand);CHKERRQ(ierr);
370     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
371     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
372     ierr = MatMult(J, du, df);CHKERRQ(ierr);
373     /* Evaluate residual at u, F(u), save in vector r */
374     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
375     ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr);
376     /* Look at the convergence of our Taylor approximation as we approach u */
377     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
378     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
379     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
380     ierr = VecDuplicate(u, &uhat_t);CHKERRQ(ierr);
381     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
382     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
383       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
384       ierr = VecWAXPY(uhat_t, h*shift, du, u_t);CHKERRQ(ierr);
385       /* F(\hat u, \hat u_t) \approx F(u, u_t) + J(u, u_t) (uhat - u) + J_t(u, u_t) (uhat_t - u_t) = F(u) + h * J(u) du + h * shift * J_t(u) du = F(u) + h F' du */
386       ierr = TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE);CHKERRQ(ierr);
387       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
388       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
389 
390       es[Nv] = PetscLog10Real(errors[Nv]);
391       hs[Nv] = PetscLog10Real(h);
392     }
393     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
394     ierr = VecDestroy(&uhat_t);CHKERRQ(ierr);
395     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
396     ierr = VecDestroy(&df);CHKERRQ(ierr);
397     ierr = VecDestroy(&r);CHKERRQ(ierr);
398     ierr = VecDestroy(&du);CHKERRQ(ierr);
399     for (v = 0; v < Nv; ++v) {
400       if ((tol >= 0) && (errors[v] > tol)) break;
401       else if (errors[v] > PETSC_SMALL)    break;
402     }
403     if (v == Nv) isLin = PETSC_TRUE;
404     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
405     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
406     /* Slope should be about 2 */
407     if (tol >= 0) {
408       if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
409     } else if (isLinear || convRate) {
410       if (isLinear) *isLinear = isLin;
411       if (convRate) *convRate = slope;
412     } else {
413       if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
414       else        {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
415     }
416   }
417   ierr = MatDestroy(&J);CHKERRQ(ierr);
418   PetscFunctionReturn(0);
419 }
420 
421 /*@C
422   DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
423 
424   Input Parameters:
425 + ts - the TS object
426 - u  - representative TS vector
427 
428   Note: The user must call PetscDSSetExactSolution() beforehand
429 
430   Level: developer
431 @*/
432 PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u)
433 {
434   DM             dm;
435   SNES           snes;
436   Vec            sol, u_t;
437   PetscReal      t;
438   PetscBool      check;
439   PetscErrorCode ierr;
440 
441   PetscFunctionBegin;
442   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
443   if (!check) PetscFunctionReturn(0);
444   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
445   ierr = VecCopy(u, sol);CHKERRQ(ierr);
446   ierr = TSSetSolution(ts, u);CHKERRQ(ierr);
447   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
448   ierr = TSSetUp(ts);CHKERRQ(ierr);
449   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
450   ierr = SNESSetSolution(snes, u);CHKERRQ(ierr);
451 
452   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
453   ierr = DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL);CHKERRQ(ierr);
454   ierr = DMGetGlobalVector(dm, &u_t);CHKERRQ(ierr);
455   ierr = DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL);CHKERRQ(ierr);
456   ierr = DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL);CHKERRQ(ierr);
457   ierr = DMRestoreGlobalVector(dm, &u_t);CHKERRQ(ierr);
458 
459   ierr = VecDestroy(&sol);CHKERRQ(ierr);
460   PetscFunctionReturn(0);
461 }
462