xref: /petsc/src/ts/interface/ts.c (revision da4996e5b0f5ea5bcee36c2408081f808bf49bc0)
1 #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
2 #include <petscdmda.h>
3 #include <petscdmshell.h>
4 #include <petscdmplex.h>  // For TSSetFromOptions()
5 #include <petscdmswarm.h> // For TSSetFromOptions()
6 #include <petscviewer.h>
7 #include <petscdraw.h>
8 #include <petscconvest.h>
9 
10 /* Logging support */
11 PetscClassId  TS_CLASSID, DMTS_CLASSID;
12 PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
13 
14 const char *const TSExactFinalTimeOptions[] = {"UNSPECIFIED", "STEPOVER", "INTERPOLATE", "MATCHSTEP", "TSExactFinalTimeOption", "TS_EXACTFINALTIME_", NULL};
15 
16 static PetscErrorCode TSAdaptSetDefaultType(TSAdapt adapt, TSAdaptType default_type)
17 {
18   PetscFunctionBegin;
19   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
20   PetscAssertPointer(default_type, 2);
21   if (!((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, default_type));
22   PetscFunctionReturn(PETSC_SUCCESS);
23 }
24 
25 /*@
26   TSSetFromOptions - Sets various `TS` parameters from the options database
27 
28   Collective
29 
30   Input Parameter:
31 . ts - the `TS` context obtained from `TSCreate()`
32 
33   Options Database Keys:
34 + -ts_type <type>                                                    - EULER, BEULER, SUNDIALS, PSEUDO, CN, RK, THETA, ALPHA, GLLE,  SSP, GLEE, BSYMP, IRK, see `TSType`
35 . -ts_save_trajectory                                                - checkpoint the solution at each time-step
36 . -ts_max_time <time>                                                - maximum time to compute to
37 . -ts_time_span <t0,...tf>                                           - sets the time span, solutions are computed and stored for each indicated time, init_time and max_time are set
38 . -ts_eval_times <t0,...tn>                                          - time points where solutions are computed and stored for each indicated time
39 . -ts_max_steps <steps>                                              - maximum number of time-steps to take
40 . -ts_init_time <time>                                               - initial time to start computation
41 . -ts_final_time <time>                                              - final time to compute to (deprecated: use `-ts_max_time`)
42 . -ts_dt <dt>                                                        - initial time step
43 . -ts_exact_final_time <stepover,interpolate,matchstep>              - whether to stop at the exact given final time and how to compute the solution at that time
44 . -ts_max_snes_failures <maxfailures>                                - Maximum number of nonlinear solve failures allowed
45 . -ts_max_reject <maxrejects>                                        - Maximum number of step rejections before step fails
46 . -ts_error_if_step_fails <true,false>                               - Error if no step succeeds
47 . -ts_rtol <rtol>                                                    - relative tolerance for local truncation error
48 . -ts_atol <atol>                                                    - Absolute tolerance for local truncation error
49 . -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view               - test the Jacobian at each iteration against finite difference with RHS function
50 . -ts_rhs_jacobian_test_mult_transpose                               - test the Jacobian at each iteration against finite difference with RHS function
51 . -ts_adjoint_solve <yes,no>                                         - After solving the ODE/DAE solve the adjoint problem (requires `-ts_save_trajectory`)
52 . -ts_fd_color                                                       - Use finite differences with coloring to compute IJacobian
53 . -ts_monitor                                                        - print information at each timestep
54 . -ts_monitor_cancel                                                 - Cancel all monitors
55 . -ts_monitor_lg_solution                                            - Monitor solution graphically
56 . -ts_monitor_lg_error                                               - Monitor error graphically
57 . -ts_monitor_error                                                  - Monitors norm of error
58 . -ts_monitor_lg_timestep                                            - Monitor timestep size graphically
59 . -ts_monitor_lg_timestep_log                                        - Monitor log timestep size graphically
60 . -ts_monitor_lg_snes_iterations                                     - Monitor number nonlinear iterations for each timestep graphically
61 . -ts_monitor_lg_ksp_iterations                                      - Monitor number nonlinear iterations for each timestep graphically
62 . -ts_monitor_sp_eig                                                 - Monitor eigenvalues of linearized operator graphically
63 . -ts_monitor_draw_solution                                          - Monitor solution graphically
64 . -ts_monitor_draw_solution_phase  <xleft,yleft,xright,yright>       - Monitor solution graphically with phase diagram, requires problem with exactly 2 degrees of freedom
65 . -ts_monitor_draw_error                                             - Monitor error graphically, requires use to have provided TSSetSolutionFunction()
66 . -ts_monitor_solution [ascii binary draw][:filename][:viewerformat] - monitors the solution at each timestep
67 . -ts_monitor_solution_interval <interval>                           - output once every interval (default=1) time steps. Use -1 to only output at the end of the simulation
68 . -ts_monitor_solution_vtk <filename.vts,filename.vtu>               - Save each time step to a binary file, use filename-%%03" PetscInt_FMT ".vts (filename-%%03" PetscInt_FMT ".vtu)
69 . -ts_monitor_solution_vtk_interval <interval>                       - output once every interval (default=1) time steps. Use -1 to only output at the end of the simulation
70 - -ts_monitor_envelope                                               - determine maximum and minimum value of each component of the solution over the solution time
71 
72   Level: beginner
73 
74   Notes:
75   See `SNESSetFromOptions()` and `KSPSetFromOptions()` for how to control the nonlinear and linear solves used by the time-stepper.
76 
77   Certain `SNES` options get reset for each new nonlinear solver, for example `-snes_lag_jacobian its` and `-snes_lag_preconditioner its`, in order
78   to retain them over the multiple nonlinear solves that `TS` uses you must also provide `-snes_lag_jacobian_persists true` and
79   `-snes_lag_preconditioner_persists true`
80 
81   Developer Notes:
82   We should unify all the -ts_monitor options in the way that -xxx_view has been unified
83 
84 .seealso: [](ch_ts), `TS`, `TSGetType()`
85 @*/
86 PetscErrorCode TSSetFromOptions(TS ts)
87 {
88   PetscBool              opt, flg, tflg;
89   char                   monfilename[PETSC_MAX_PATH_LEN];
90   PetscReal              time_step, eval_times[100];
91   PetscInt               num_eval_times = PETSC_STATIC_ARRAY_LENGTH(eval_times);
92   TSExactFinalTimeOption eftopt;
93   char                   dir[16];
94   TSIFunctionFn         *ifun;
95   const char            *defaultType;
96   char                   typeName[256];
97 
98   PetscFunctionBegin;
99   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
100 
101   PetscCall(TSRegisterAll());
102   PetscCall(TSGetIFunction(ts, NULL, &ifun, NULL));
103 
104   PetscObjectOptionsBegin((PetscObject)ts);
105   if (((PetscObject)ts)->type_name) defaultType = ((PetscObject)ts)->type_name;
106   else defaultType = ifun ? TSBEULER : TSEULER;
107   PetscCall(PetscOptionsFList("-ts_type", "TS method", "TSSetType", TSList, defaultType, typeName, 256, &opt));
108   if (opt) PetscCall(TSSetType(ts, typeName));
109   else PetscCall(TSSetType(ts, defaultType));
110 
111   /* Handle generic TS options */
112   PetscCall(PetscOptionsDeprecated("-ts_final_time", "-ts_max_time", "3.10", NULL));
113   PetscCall(PetscOptionsReal("-ts_max_time", "Maximum time to run to", "TSSetMaxTime", ts->max_time, &ts->max_time, NULL));
114   PetscCall(PetscOptionsRealArray("-ts_time_span", "Time span", "TSSetTimeSpan", eval_times, &num_eval_times, &flg));
115   if (flg) PetscCall(TSSetTimeSpan(ts, num_eval_times, eval_times));
116   num_eval_times = PETSC_STATIC_ARRAY_LENGTH(eval_times);
117   PetscCall(PetscOptionsRealArray("-ts_eval_times", "Evaluation time points", "TSSetEvaluationTimes", eval_times, &num_eval_times, &opt));
118   PetscCheck(flg != opt || (!flg && !opt), PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "May not provide -ts_time_span and -ts_eval_times simultaneously");
119   if (opt) PetscCall(TSSetEvaluationTimes(ts, num_eval_times, eval_times));
120   PetscCall(PetscOptionsInt("-ts_max_steps", "Maximum number of time steps", "TSSetMaxSteps", ts->max_steps, &ts->max_steps, NULL));
121   PetscCall(PetscOptionsReal("-ts_init_time", "Initial time", "TSSetTime", ts->ptime, &ts->ptime, NULL));
122   PetscCall(PetscOptionsReal("-ts_dt", "Initial time step", "TSSetTimeStep", ts->time_step, &time_step, &flg));
123   if (flg) PetscCall(TSSetTimeStep(ts, time_step));
124   PetscCall(PetscOptionsEnum("-ts_exact_final_time", "Option for handling of final time step", "TSSetExactFinalTime", TSExactFinalTimeOptions, (PetscEnum)ts->exact_final_time, (PetscEnum *)&eftopt, &flg));
125   if (flg) PetscCall(TSSetExactFinalTime(ts, eftopt));
126   PetscCall(PetscOptionsInt("-ts_max_snes_failures", "Maximum number of nonlinear solve failures", "TSSetMaxSNESFailures", ts->max_snes_failures, &ts->max_snes_failures, &flg));
127   if (flg) PetscCall(TSSetMaxSNESFailures(ts, ts->max_snes_failures));
128   PetscCall(PetscOptionsInt("-ts_max_reject", "Maximum number of step rejections before step fails", "TSSetMaxStepRejections", ts->max_reject, &ts->max_reject, &flg));
129   if (flg) PetscCall(TSSetMaxStepRejections(ts, ts->max_reject));
130   PetscCall(PetscOptionsBool("-ts_error_if_step_fails", "Error if no step succeeds", "TSSetErrorIfStepFails", ts->errorifstepfailed, &ts->errorifstepfailed, NULL));
131   PetscCall(PetscOptionsBoundedReal("-ts_rtol", "Relative tolerance for local truncation error", "TSSetTolerances", ts->rtol, &ts->rtol, NULL, 0));
132   PetscCall(PetscOptionsBoundedReal("-ts_atol", "Absolute tolerance for local truncation error", "TSSetTolerances", ts->atol, &ts->atol, NULL, 0));
133 
134   PetscCall(PetscOptionsBool("-ts_rhs_jacobian_test_mult", "Test the RHS Jacobian for consistency with RHS at each solve ", "None", ts->testjacobian, &ts->testjacobian, NULL));
135   PetscCall(PetscOptionsBool("-ts_rhs_jacobian_test_mult_transpose", "Test the RHS Jacobian transpose for consistency with RHS at each solve ", "None", ts->testjacobiantranspose, &ts->testjacobiantranspose, NULL));
136   PetscCall(PetscOptionsBool("-ts_use_splitrhsfunction", "Use the split RHS function for multirate solvers ", "TSSetUseSplitRHSFunction", ts->use_splitrhsfunction, &ts->use_splitrhsfunction, NULL));
137 #if defined(PETSC_HAVE_SAWS)
138   {
139     PetscBool set;
140     flg = PETSC_FALSE;
141     PetscCall(PetscOptionsBool("-ts_saws_block", "Block for SAWs memory snooper at end of TSSolve", "PetscObjectSAWsBlock", ((PetscObject)ts)->amspublishblock, &flg, &set));
142     if (set) PetscCall(PetscObjectSAWsSetBlock((PetscObject)ts, flg));
143   }
144 #endif
145 
146   /* Monitor options */
147   PetscCall(PetscOptionsInt("-ts_monitor_frequency", "Number of time steps between monitor output", "TSMonitorSetFrequency", ts->monitorFrequency, &ts->monitorFrequency, NULL));
148   PetscCall(TSMonitorSetFromOptions(ts, "-ts_monitor", "Monitor time and timestep size", "TSMonitorDefault", TSMonitorDefault, NULL));
149   PetscCall(TSMonitorSetFromOptions(ts, "-ts_monitor_extreme", "Monitor extreme values of the solution", "TSMonitorExtreme", TSMonitorExtreme, NULL));
150   PetscCall(TSMonitorSetFromOptions(ts, "-ts_monitor_solution", "View the solution at each timestep", "TSMonitorSolution", TSMonitorSolution, NULL));
151   PetscCall(TSMonitorSetFromOptions(ts, "-ts_dmswarm_monitor_moments", "Monitor moments of particle distribution", "TSDMSwarmMonitorMoments", TSDMSwarmMonitorMoments, NULL));
152 
153   PetscCall(PetscOptionsString("-ts_monitor_python", "Use Python function", "TSMonitorSet", NULL, monfilename, sizeof(monfilename), &flg));
154   if (flg) PetscCall(PetscPythonMonitorSet((PetscObject)ts, monfilename));
155 
156   PetscCall(PetscOptionsName("-ts_monitor_lg_solution", "Monitor solution graphically", "TSMonitorLGSolution", &opt));
157   if (opt) {
158     PetscInt  howoften = 1;
159     DM        dm;
160     PetscBool net;
161 
162     PetscCall(PetscOptionsInt("-ts_monitor_lg_solution", "Monitor solution graphically", "TSMonitorLGSolution", howoften, &howoften, NULL));
163     PetscCall(TSGetDM(ts, &dm));
164     PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMNETWORK, &net));
165     if (net) {
166       TSMonitorLGCtxNetwork ctx;
167       PetscCall(TSMonitorLGCtxNetworkCreate(ts, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 600, 400, howoften, &ctx));
168       PetscCall(TSMonitorSet(ts, TSMonitorLGCtxNetworkSolution, ctx, (PetscCtxDestroyFn *)TSMonitorLGCtxNetworkDestroy));
169       PetscCall(PetscOptionsBool("-ts_monitor_lg_solution_semilogy", "Plot the solution with a semi-log axis", "", ctx->semilogy, &ctx->semilogy, NULL));
170     } else {
171       TSMonitorLGCtx ctx;
172       PetscCall(TSMonitorLGCtxCreate(PETSC_COMM_SELF, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, howoften, &ctx));
173       PetscCall(TSMonitorSet(ts, TSMonitorLGSolution, ctx, (PetscCtxDestroyFn *)TSMonitorLGCtxDestroy));
174     }
175   }
176 
177   PetscCall(PetscOptionsName("-ts_monitor_lg_error", "Monitor error graphically", "TSMonitorLGError", &opt));
178   if (opt) {
179     TSMonitorLGCtx ctx;
180     PetscInt       howoften = 1;
181 
182     PetscCall(PetscOptionsInt("-ts_monitor_lg_error", "Monitor error graphically", "TSMonitorLGError", howoften, &howoften, NULL));
183     PetscCall(TSMonitorLGCtxCreate(PETSC_COMM_SELF, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, howoften, &ctx));
184     PetscCall(TSMonitorSet(ts, TSMonitorLGError, ctx, (PetscCtxDestroyFn *)TSMonitorLGCtxDestroy));
185   }
186   PetscCall(TSMonitorSetFromOptions(ts, "-ts_monitor_error", "View the error at each timestep", "TSMonitorError", TSMonitorError, NULL));
187 
188   PetscCall(PetscOptionsName("-ts_monitor_lg_timestep", "Monitor timestep size graphically", "TSMonitorLGTimeStep", &opt));
189   if (opt) {
190     TSMonitorLGCtx ctx;
191     PetscInt       howoften = 1;
192 
193     PetscCall(PetscOptionsInt("-ts_monitor_lg_timestep", "Monitor timestep size graphically", "TSMonitorLGTimeStep", howoften, &howoften, NULL));
194     PetscCall(TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, howoften, &ctx));
195     PetscCall(TSMonitorSet(ts, TSMonitorLGTimeStep, ctx, (PetscCtxDestroyFn *)TSMonitorLGCtxDestroy));
196   }
197   PetscCall(PetscOptionsName("-ts_monitor_lg_timestep_log", "Monitor log timestep size graphically", "TSMonitorLGTimeStep", &opt));
198   if (opt) {
199     TSMonitorLGCtx ctx;
200     PetscInt       howoften = 1;
201 
202     PetscCall(PetscOptionsInt("-ts_monitor_lg_timestep_log", "Monitor log timestep size graphically", "TSMonitorLGTimeStep", howoften, &howoften, NULL));
203     PetscCall(TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, howoften, &ctx));
204     PetscCall(TSMonitorSet(ts, TSMonitorLGTimeStep, ctx, (PetscCtxDestroyFn *)TSMonitorLGCtxDestroy));
205     ctx->semilogy = PETSC_TRUE;
206   }
207 
208   PetscCall(PetscOptionsName("-ts_monitor_lg_snes_iterations", "Monitor number nonlinear iterations for each timestep graphically", "TSMonitorLGSNESIterations", &opt));
209   if (opt) {
210     TSMonitorLGCtx ctx;
211     PetscInt       howoften = 1;
212 
213     PetscCall(PetscOptionsInt("-ts_monitor_lg_snes_iterations", "Monitor number nonlinear iterations for each timestep graphically", "TSMonitorLGSNESIterations", howoften, &howoften, NULL));
214     PetscCall(TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, howoften, &ctx));
215     PetscCall(TSMonitorSet(ts, TSMonitorLGSNESIterations, ctx, (PetscCtxDestroyFn *)TSMonitorLGCtxDestroy));
216   }
217   PetscCall(PetscOptionsName("-ts_monitor_lg_ksp_iterations", "Monitor number nonlinear iterations for each timestep graphically", "TSMonitorLGKSPIterations", &opt));
218   if (opt) {
219     TSMonitorLGCtx ctx;
220     PetscInt       howoften = 1;
221 
222     PetscCall(PetscOptionsInt("-ts_monitor_lg_ksp_iterations", "Monitor number nonlinear iterations for each timestep graphically", "TSMonitorLGKSPIterations", howoften, &howoften, NULL));
223     PetscCall(TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, howoften, &ctx));
224     PetscCall(TSMonitorSet(ts, TSMonitorLGKSPIterations, ctx, (PetscCtxDestroyFn *)TSMonitorLGCtxDestroy));
225   }
226   PetscCall(PetscOptionsName("-ts_monitor_sp_eig", "Monitor eigenvalues of linearized operator graphically", "TSMonitorSPEig", &opt));
227   if (opt) {
228     TSMonitorSPEigCtx ctx;
229     PetscInt          howoften = 1;
230 
231     PetscCall(PetscOptionsInt("-ts_monitor_sp_eig", "Monitor eigenvalues of linearized operator graphically", "TSMonitorSPEig", howoften, &howoften, NULL));
232     PetscCall(TSMonitorSPEigCtxCreate(PETSC_COMM_SELF, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx));
233     PetscCall(TSMonitorSet(ts, TSMonitorSPEig, ctx, (PetscCtxDestroyFn *)TSMonitorSPEigCtxDestroy));
234   }
235   PetscCall(PetscOptionsName("-ts_monitor_sp_swarm", "Display particle phase space from the DMSwarm", "TSMonitorSPSwarm", &opt));
236   if (opt) {
237     TSMonitorSPCtx ctx;
238     PetscInt       howoften = 1, retain = 0;
239     PetscBool      phase = PETSC_TRUE, create = PETSC_TRUE, multispecies = PETSC_FALSE;
240 
241     for (PetscInt i = 0; i < ts->numbermonitors; ++i)
242       if (ts->monitor[i] == TSMonitorSPSwarmSolution) {
243         create = PETSC_FALSE;
244         break;
245       }
246     if (create) {
247       PetscCall(PetscOptionsInt("-ts_monitor_sp_swarm", "Display particles phase space from the DMSwarm", "TSMonitorSPSwarm", howoften, &howoften, NULL));
248       PetscCall(PetscOptionsInt("-ts_monitor_sp_swarm_retain", "Retain n points plotted to show trajectory, -1 for all points", "TSMonitorSPSwarm", retain, &retain, NULL));
249       PetscCall(PetscOptionsBool("-ts_monitor_sp_swarm_phase", "Plot in phase space rather than coordinate space", "TSMonitorSPSwarm", phase, &phase, NULL));
250       PetscCall(PetscOptionsBool("-ts_monitor_sp_swarm_multi_species", "Color particles by particle species", "TSMonitorSPSwarm", multispecies, &multispecies, NULL));
251       PetscCall(TSMonitorSPCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, retain, phase, multispecies, &ctx));
252       PetscCall(TSMonitorSet(ts, TSMonitorSPSwarmSolution, ctx, (PetscCtxDestroyFn *)TSMonitorSPCtxDestroy));
253     }
254   }
255   PetscCall(PetscOptionsName("-ts_monitor_hg_swarm", "Display particle histogram from the DMSwarm", "TSMonitorHGSwarm", &opt));
256   if (opt) {
257     TSMonitorHGCtx ctx;
258     PetscInt       howoften = 1, Ns = 1;
259     PetscBool      velocity = PETSC_FALSE, create = PETSC_TRUE;
260 
261     for (PetscInt i = 0; i < ts->numbermonitors; ++i)
262       if (ts->monitor[i] == TSMonitorHGSwarmSolution) {
263         create = PETSC_FALSE;
264         break;
265       }
266     if (create) {
267       DM       sw, dm;
268       PetscInt Nc, Nb;
269 
270       PetscCall(TSGetDM(ts, &sw));
271       PetscCall(DMSwarmGetCellDM(sw, &dm));
272       PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &Nc));
273       Nb = PetscMin(20, PetscMax(10, Nc));
274       PetscCall(PetscOptionsInt("-ts_monitor_hg_swarm", "Display particles histogram from the DMSwarm", "TSMonitorHGSwarm", howoften, &howoften, NULL));
275       PetscCall(PetscOptionsBool("-ts_monitor_hg_swarm_velocity", "Plot in velocity space rather than coordinate space", "TSMonitorHGSwarm", velocity, &velocity, NULL));
276       PetscCall(PetscOptionsInt("-ts_monitor_hg_swarm_species", "Number of species to histogram", "TSMonitorHGSwarm", Ns, &Ns, NULL));
277       PetscCall(PetscOptionsInt("-ts_monitor_hg_swarm_bins", "Number of histogram bins", "TSMonitorHGSwarm", Nb, &Nb, NULL));
278       PetscCall(TSMonitorHGCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, Ns, Nb, velocity, &ctx));
279       PetscCall(TSMonitorSet(ts, TSMonitorHGSwarmSolution, ctx, (PetscCtxDestroyFn *)TSMonitorHGCtxDestroy));
280     }
281   }
282   opt = PETSC_FALSE;
283   PetscCall(PetscOptionsName("-ts_monitor_draw_solution", "Monitor solution graphically", "TSMonitorDrawSolution", &opt));
284   if (opt) {
285     TSMonitorDrawCtx ctx;
286     PetscInt         howoften = 1;
287 
288     PetscCall(PetscOptionsInt("-ts_monitor_draw_solution", "Monitor solution graphically", "TSMonitorDrawSolution", howoften, &howoften, NULL));
289     PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, "Computed Solution", PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx));
290     PetscCall(TSMonitorSet(ts, TSMonitorDrawSolution, ctx, (PetscCtxDestroyFn *)TSMonitorDrawCtxDestroy));
291   }
292   opt = PETSC_FALSE;
293   PetscCall(PetscOptionsName("-ts_monitor_draw_solution_phase", "Monitor solution graphically", "TSMonitorDrawSolutionPhase", &opt));
294   if (opt) {
295     TSMonitorDrawCtx ctx;
296     PetscReal        bounds[4];
297     PetscInt         n = 4;
298     PetscDraw        draw;
299     PetscDrawAxis    axis;
300 
301     PetscCall(PetscOptionsRealArray("-ts_monitor_draw_solution_phase", "Monitor solution graphically", "TSMonitorDrawSolutionPhase", bounds, &n, NULL));
302     PetscCheck(n == 4, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Must provide bounding box of phase field");
303     PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, 1, &ctx));
304     PetscCall(PetscViewerDrawGetDraw(ctx->viewer, 0, &draw));
305     PetscCall(PetscViewerDrawGetDrawAxis(ctx->viewer, 0, &axis));
306     PetscCall(PetscDrawAxisSetLimits(axis, bounds[0], bounds[2], bounds[1], bounds[3]));
307     PetscCall(PetscDrawAxisSetLabels(axis, "Phase Diagram", "Variable 1", "Variable 2"));
308     PetscCall(TSMonitorSet(ts, TSMonitorDrawSolutionPhase, ctx, (PetscCtxDestroyFn *)TSMonitorDrawCtxDestroy));
309   }
310   opt = PETSC_FALSE;
311   PetscCall(PetscOptionsName("-ts_monitor_draw_error", "Monitor error graphically", "TSMonitorDrawError", &opt));
312   if (opt) {
313     TSMonitorDrawCtx ctx;
314     PetscInt         howoften = 1;
315 
316     PetscCall(PetscOptionsInt("-ts_monitor_draw_error", "Monitor error graphically", "TSMonitorDrawError", howoften, &howoften, NULL));
317     PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, "Error", PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx));
318     PetscCall(TSMonitorSet(ts, TSMonitorDrawError, ctx, (PetscCtxDestroyFn *)TSMonitorDrawCtxDestroy));
319   }
320   opt = PETSC_FALSE;
321   PetscCall(PetscOptionsName("-ts_monitor_draw_solution_function", "Monitor solution provided by TSMonitorSetSolutionFunction() graphically", "TSMonitorDrawSolutionFunction", &opt));
322   if (opt) {
323     TSMonitorDrawCtx ctx;
324     PetscInt         howoften = 1;
325 
326     PetscCall(PetscOptionsInt("-ts_monitor_draw_solution_function", "Monitor solution provided by TSMonitorSetSolutionFunction() graphically", "TSMonitorDrawSolutionFunction", howoften, &howoften, NULL));
327     PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, "Solution provided by user function", PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx));
328     PetscCall(TSMonitorSet(ts, TSMonitorDrawSolutionFunction, ctx, (PetscCtxDestroyFn *)TSMonitorDrawCtxDestroy));
329   }
330 
331   opt = PETSC_FALSE;
332   PetscCall(PetscOptionsString("-ts_monitor_solution_vtk", "Save each time step to a binary file, use filename-%%03" PetscInt_FMT ".vts", "TSMonitorSolutionVTK", NULL, monfilename, sizeof(monfilename), &flg));
333   if (flg) {
334     TSMonitorVTKCtx ctx;
335 
336     PetscCall(TSMonitorSolutionVTKCtxCreate(monfilename, &ctx));
337     PetscCall(PetscOptionsInt("-ts_monitor_solution_vtk_interval", "Save every interval time step (-1 for last step only)", NULL, ctx->interval, &ctx->interval, NULL));
338     PetscCall(TSMonitorSet(ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))TSMonitorSolutionVTK, ctx, (PetscCtxDestroyFn *)TSMonitorSolutionVTKDestroy));
339   }
340 
341   PetscCall(PetscOptionsString("-ts_monitor_dmda_ray", "Display a ray of the solution", "None", "y=0", dir, sizeof(dir), &flg));
342   if (flg) {
343     TSMonitorDMDARayCtx *rayctx;
344     int                  ray = 0;
345     DMDirection          ddir;
346     DM                   da;
347     PetscMPIInt          rank;
348 
349     PetscCheck(dir[1] == '=', PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Unknown ray %s", dir);
350     if (dir[0] == 'x') ddir = DM_X;
351     else if (dir[0] == 'y') ddir = DM_Y;
352     else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Unknown ray %s", dir);
353     sscanf(dir + 2, "%d", &ray);
354 
355     PetscCall(PetscInfo(ts, "Displaying DMDA ray %c = %d\n", dir[0], ray));
356     PetscCall(PetscNew(&rayctx));
357     PetscCall(TSGetDM(ts, &da));
358     PetscCall(DMDAGetRay(da, ddir, ray, &rayctx->ray, &rayctx->scatter));
359     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ts), &rank));
360     if (rank == 0) PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, NULL, NULL, 0, 0, 600, 300, &rayctx->viewer));
361     rayctx->lgctx = NULL;
362     PetscCall(TSMonitorSet(ts, TSMonitorDMDARay, rayctx, TSMonitorDMDARayDestroy));
363   }
364   PetscCall(PetscOptionsString("-ts_monitor_lg_dmda_ray", "Display a ray of the solution", "None", "x=0", dir, sizeof(dir), &flg));
365   if (flg) {
366     TSMonitorDMDARayCtx *rayctx;
367     int                  ray = 0;
368     DMDirection          ddir;
369     DM                   da;
370     PetscInt             howoften = 1;
371 
372     PetscCheck(dir[1] == '=', PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Malformed ray %s", dir);
373     if (dir[0] == 'x') ddir = DM_X;
374     else if (dir[0] == 'y') ddir = DM_Y;
375     else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Unknown ray direction %s", dir);
376     sscanf(dir + 2, "%d", &ray);
377 
378     PetscCall(PetscInfo(ts, "Displaying LG DMDA ray %c = %d\n", dir[0], ray));
379     PetscCall(PetscNew(&rayctx));
380     PetscCall(TSGetDM(ts, &da));
381     PetscCall(DMDAGetRay(da, ddir, ray, &rayctx->ray, &rayctx->scatter));
382     PetscCall(TSMonitorLGCtxCreate(PETSC_COMM_SELF, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 600, 400, howoften, &rayctx->lgctx));
383     PetscCall(TSMonitorSet(ts, TSMonitorLGDMDARay, rayctx, TSMonitorDMDARayDestroy));
384   }
385 
386   PetscCall(PetscOptionsName("-ts_monitor_envelope", "Monitor maximum and minimum value of each component of the solution", "TSMonitorEnvelope", &opt));
387   if (opt) {
388     TSMonitorEnvelopeCtx ctx;
389 
390     PetscCall(TSMonitorEnvelopeCtxCreate(ts, &ctx));
391     PetscCall(TSMonitorSet(ts, TSMonitorEnvelope, ctx, (PetscCtxDestroyFn *)TSMonitorEnvelopeCtxDestroy));
392   }
393   flg = PETSC_FALSE;
394   PetscCall(PetscOptionsBool("-ts_monitor_cancel", "Remove all monitors", "TSMonitorCancel", flg, &flg, &opt));
395   if (opt && flg) PetscCall(TSMonitorCancel(ts));
396 
397   flg = PETSC_FALSE;
398   PetscCall(PetscOptionsBool("-ts_fd_color", "Use finite differences with coloring to compute IJacobian", "TSComputeIJacobianDefaultColor", flg, &flg, NULL));
399   if (flg) {
400     DM dm;
401 
402     PetscCall(TSGetDM(ts, &dm));
403     PetscCall(DMTSUnsetIJacobianContext_Internal(dm));
404     PetscCall(TSSetIJacobian(ts, NULL, NULL, TSComputeIJacobianDefaultColor, NULL));
405     PetscCall(PetscInfo(ts, "Setting default finite difference coloring Jacobian matrix\n"));
406   }
407 
408   /* Handle specific TS options */
409   PetscTryTypeMethod(ts, setfromoptions, PetscOptionsObject);
410 
411   /* Handle TSAdapt options */
412   PetscCall(TSGetAdapt(ts, &ts->adapt));
413   PetscCall(TSAdaptSetDefaultType(ts->adapt, ts->default_adapt_type));
414   PetscCall(TSAdaptSetFromOptions(ts->adapt, PetscOptionsObject));
415 
416   /* TS trajectory must be set after TS, since it may use some TS options above */
417   tflg = ts->trajectory ? PETSC_TRUE : PETSC_FALSE;
418   PetscCall(PetscOptionsBool("-ts_save_trajectory", "Save the solution at each timestep", "TSSetSaveTrajectory", tflg, &tflg, NULL));
419   if (tflg) PetscCall(TSSetSaveTrajectory(ts));
420 
421   PetscCall(TSAdjointSetFromOptions(ts, PetscOptionsObject));
422 
423   /* process any options handlers added with PetscObjectAddOptionsHandler() */
424   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)ts, PetscOptionsObject));
425   PetscOptionsEnd();
426 
427   if (ts->trajectory) PetscCall(TSTrajectorySetFromOptions(ts->trajectory, ts));
428 
429   /* why do we have to do this here and not during TSSetUp? */
430   PetscCall(TSGetSNES(ts, &ts->snes));
431   if (ts->problem_type == TS_LINEAR) {
432     PetscCall(PetscObjectTypeCompareAny((PetscObject)ts->snes, &flg, SNESKSPONLY, SNESKSPTRANSPOSEONLY, ""));
433     if (!flg) PetscCall(SNESSetType(ts->snes, SNESKSPONLY));
434   }
435   PetscCall(SNESSetFromOptions(ts->snes));
436   PetscFunctionReturn(PETSC_SUCCESS);
437 }
438 
439 /*@
440   TSGetTrajectory - Gets the trajectory from a `TS` if it exists
441 
442   Collective
443 
444   Input Parameter:
445 . ts - the `TS` context obtained from `TSCreate()`
446 
447   Output Parameter:
448 . tr - the `TSTrajectory` object, if it exists
449 
450   Level: advanced
451 
452   Note:
453   This routine should be called after all `TS` options have been set
454 
455 .seealso: [](ch_ts), `TS`, `TSTrajectory`, `TSAdjointSolve()`, `TSTrajectoryCreate()`
456 @*/
457 PetscErrorCode TSGetTrajectory(TS ts, TSTrajectory *tr)
458 {
459   PetscFunctionBegin;
460   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
461   *tr = ts->trajectory;
462   PetscFunctionReturn(PETSC_SUCCESS);
463 }
464 
465 /*@
466   TSSetSaveTrajectory - Causes the `TS` to save its solutions as it iterates forward in time in a `TSTrajectory` object
467 
468   Collective
469 
470   Input Parameter:
471 . ts - the `TS` context obtained from `TSCreate()`
472 
473   Options Database Keys:
474 + -ts_save_trajectory      - saves the trajectory to a file
475 - -ts_trajectory_type type - set trajectory type
476 
477   Level: intermediate
478 
479   Notes:
480   This routine should be called after all `TS` options have been set
481 
482   The `TSTRAJECTORYVISUALIZATION` files can be loaded into Python with $PETSC_DIR/lib/petsc/bin/PetscBinaryIOTrajectory.py and
483   MATLAB with $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m
484 
485 .seealso: [](ch_ts), `TS`, `TSTrajectory`, `TSGetTrajectory()`, `TSAdjointSolve()`
486 @*/
487 PetscErrorCode TSSetSaveTrajectory(TS ts)
488 {
489   PetscFunctionBegin;
490   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
491   if (!ts->trajectory) PetscCall(TSTrajectoryCreate(PetscObjectComm((PetscObject)ts), &ts->trajectory));
492   PetscFunctionReturn(PETSC_SUCCESS);
493 }
494 
495 /*@
496   TSResetTrajectory - Destroys and recreates the internal `TSTrajectory` object
497 
498   Collective
499 
500   Input Parameter:
501 . ts - the `TS` context obtained from `TSCreate()`
502 
503   Level: intermediate
504 
505 .seealso: [](ch_ts), `TSTrajectory`, `TSGetTrajectory()`, `TSAdjointSolve()`, `TSRemoveTrajectory()`
506 @*/
507 PetscErrorCode TSResetTrajectory(TS ts)
508 {
509   PetscFunctionBegin;
510   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
511   if (ts->trajectory) {
512     PetscCall(TSTrajectoryDestroy(&ts->trajectory));
513     PetscCall(TSTrajectoryCreate(PetscObjectComm((PetscObject)ts), &ts->trajectory));
514   }
515   PetscFunctionReturn(PETSC_SUCCESS);
516 }
517 
518 /*@
519   TSRemoveTrajectory - Destroys and removes the internal `TSTrajectory` object from a `TS`
520 
521   Collective
522 
523   Input Parameter:
524 . ts - the `TS` context obtained from `TSCreate()`
525 
526   Level: intermediate
527 
528 .seealso: [](ch_ts), `TSTrajectory`, `TSResetTrajectory()`, `TSAdjointSolve()`
529 @*/
530 PetscErrorCode TSRemoveTrajectory(TS ts)
531 {
532   PetscFunctionBegin;
533   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
534   if (ts->trajectory) PetscCall(TSTrajectoryDestroy(&ts->trajectory));
535   PetscFunctionReturn(PETSC_SUCCESS);
536 }
537 
538 /*@
539   TSComputeRHSJacobian - Computes the Jacobian matrix that has been
540   set with `TSSetRHSJacobian()`.
541 
542   Collective
543 
544   Input Parameters:
545 + ts - the `TS` context
546 . t  - current timestep
547 - U  - input vector
548 
549   Output Parameters:
550 + A - Jacobian matrix
551 - B - optional preconditioning matrix
552 
553   Level: developer
554 
555   Note:
556   Most users should not need to explicitly call this routine, as it
557   is used internally within the nonlinear solvers.
558 
559 .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `KSPSetOperators()`
560 @*/
561 PetscErrorCode TSComputeRHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B)
562 {
563   PetscObjectState Ustate;
564   PetscObjectId    Uid;
565   DM               dm;
566   DMTS             tsdm;
567   TSRHSJacobianFn *rhsjacobianfunc;
568   void            *ctx;
569   TSRHSFunctionFn *rhsfunction;
570 
571   PetscFunctionBegin;
572   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
573   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
574   PetscCheckSameComm(ts, 1, U, 3);
575   PetscCall(TSGetDM(ts, &dm));
576   PetscCall(DMGetDMTS(dm, &tsdm));
577   PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL));
578   PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobianfunc, &ctx));
579   PetscCall(PetscObjectStateGet((PetscObject)U, &Ustate));
580   PetscCall(PetscObjectGetId((PetscObject)U, &Uid));
581 
582   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.Xid == Uid && ts->rhsjacobian.Xstate == Ustate)) && (rhsfunction != TSComputeRHSFunctionLinear)) PetscFunctionReturn(PETSC_SUCCESS);
583 
584   PetscCheck(ts->rhsjacobian.shift == 0.0 || !ts->rhsjacobian.reuse, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Should not call TSComputeRHSJacobian() on a shifted matrix (shift=%lf) when RHSJacobian is reusable.", (double)ts->rhsjacobian.shift);
585   if (rhsjacobianfunc) {
586     PetscCall(PetscLogEventBegin(TS_JacobianEval, U, ts, A, B));
587     PetscCallBack("TS callback Jacobian", (*rhsjacobianfunc)(ts, t, U, A, B, ctx));
588     ts->rhsjacs++;
589     PetscCall(PetscLogEventEnd(TS_JacobianEval, U, ts, A, B));
590   } else {
591     PetscCall(MatZeroEntries(A));
592     if (B && A != B) PetscCall(MatZeroEntries(B));
593   }
594   ts->rhsjacobian.time  = t;
595   ts->rhsjacobian.shift = 0;
596   ts->rhsjacobian.scale = 1.;
597   PetscCall(PetscObjectGetId((PetscObject)U, &ts->rhsjacobian.Xid));
598   PetscCall(PetscObjectStateGet((PetscObject)U, &ts->rhsjacobian.Xstate));
599   PetscFunctionReturn(PETSC_SUCCESS);
600 }
601 
602 /*@
603   TSComputeRHSFunction - Evaluates the right-hand-side function for a `TS`
604 
605   Collective
606 
607   Input Parameters:
608 + ts - the `TS` context
609 . t  - current time
610 - U  - state vector
611 
612   Output Parameter:
613 . y - right-hand side
614 
615   Level: developer
616 
617   Note:
618   Most users should not need to explicitly call this routine, as it
619   is used internally within the nonlinear solvers.
620 
621 .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `TSComputeIFunction()`
622 @*/
623 PetscErrorCode TSComputeRHSFunction(TS ts, PetscReal t, Vec U, Vec y)
624 {
625   TSRHSFunctionFn *rhsfunction;
626   TSIFunctionFn   *ifunction;
627   void            *ctx;
628   DM               dm;
629 
630   PetscFunctionBegin;
631   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
632   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
633   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
634   PetscCall(TSGetDM(ts, &dm));
635   PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, &ctx));
636   PetscCall(DMTSGetIFunction(dm, &ifunction, NULL));
637 
638   PetscCheck(rhsfunction || ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSFunction() and / or TSSetIFunction()");
639 
640   if (rhsfunction) {
641     PetscCall(PetscLogEventBegin(TS_FunctionEval, U, ts, y, 0));
642     PetscCall(VecLockReadPush(U));
643     PetscCallBack("TS callback right-hand-side", (*rhsfunction)(ts, t, U, y, ctx));
644     PetscCall(VecLockReadPop(U));
645     ts->rhsfuncs++;
646     PetscCall(PetscLogEventEnd(TS_FunctionEval, U, ts, y, 0));
647   } else PetscCall(VecZeroEntries(y));
648   PetscFunctionReturn(PETSC_SUCCESS);
649 }
650 
651 /*@
652   TSComputeSolutionFunction - Evaluates the solution function.
653 
654   Collective
655 
656   Input Parameters:
657 + ts - the `TS` context
658 - t  - current time
659 
660   Output Parameter:
661 . U - the solution
662 
663   Level: developer
664 
665 .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `TSSetRHSFunction()`, `TSComputeIFunction()`
666 @*/
667 PetscErrorCode TSComputeSolutionFunction(TS ts, PetscReal t, Vec U)
668 {
669   TSSolutionFn *solutionfunction;
670   void         *ctx;
671   DM            dm;
672 
673   PetscFunctionBegin;
674   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
675   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
676   PetscCall(TSGetDM(ts, &dm));
677   PetscCall(DMTSGetSolutionFunction(dm, &solutionfunction, &ctx));
678   if (solutionfunction) PetscCallBack("TS callback solution", (*solutionfunction)(ts, t, U, ctx));
679   PetscFunctionReturn(PETSC_SUCCESS);
680 }
681 /*@
682   TSComputeForcingFunction - Evaluates the forcing function.
683 
684   Collective
685 
686   Input Parameters:
687 + ts - the `TS` context
688 - t  - current time
689 
690   Output Parameter:
691 . U - the function value
692 
693   Level: developer
694 
695 .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `TSSetRHSFunction()`, `TSComputeIFunction()`
696 @*/
697 PetscErrorCode TSComputeForcingFunction(TS ts, PetscReal t, Vec U)
698 {
699   void        *ctx;
700   DM           dm;
701   TSForcingFn *forcing;
702 
703   PetscFunctionBegin;
704   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
705   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
706   PetscCall(TSGetDM(ts, &dm));
707   PetscCall(DMTSGetForcingFunction(dm, &forcing, &ctx));
708 
709   if (forcing) PetscCallBack("TS callback forcing function", (*forcing)(ts, t, U, ctx));
710   PetscFunctionReturn(PETSC_SUCCESS);
711 }
712 
713 PetscErrorCode TSGetRHSMats_Private(TS ts, Mat *Arhs, Mat *Brhs)
714 {
715   Mat            A, B;
716   TSIJacobianFn *ijacobian;
717 
718   PetscFunctionBegin;
719   if (Arhs) *Arhs = NULL;
720   if (Brhs) *Brhs = NULL;
721   PetscCall(TSGetIJacobian(ts, &A, &B, &ijacobian, NULL));
722   if (Arhs) {
723     if (!ts->Arhs) {
724       if (ijacobian) {
725         PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &ts->Arhs));
726         PetscCall(TSSetMatStructure(ts, SAME_NONZERO_PATTERN));
727       } else {
728         ts->Arhs = A;
729         PetscCall(PetscObjectReference((PetscObject)A));
730       }
731     } else {
732       PetscBool flg;
733       PetscCall(SNESGetUseMatrixFree(ts->snes, NULL, &flg));
734       /* Handle case where user provided only RHSJacobian and used -snes_mf_operator */
735       if (flg && !ijacobian && ts->Arhs == ts->Brhs) {
736         PetscCall(PetscObjectDereference((PetscObject)ts->Arhs));
737         ts->Arhs = A;
738         PetscCall(PetscObjectReference((PetscObject)A));
739       }
740     }
741     *Arhs = ts->Arhs;
742   }
743   if (Brhs) {
744     if (!ts->Brhs) {
745       if (A != B) {
746         if (ijacobian) {
747           PetscCall(MatDuplicate(B, MAT_DO_NOT_COPY_VALUES, &ts->Brhs));
748         } else {
749           ts->Brhs = B;
750           PetscCall(PetscObjectReference((PetscObject)B));
751         }
752       } else {
753         PetscCall(PetscObjectReference((PetscObject)ts->Arhs));
754         ts->Brhs = ts->Arhs;
755       }
756     }
757     *Brhs = ts->Brhs;
758   }
759   PetscFunctionReturn(PETSC_SUCCESS);
760 }
761 
762 /*@
763   TSComputeIFunction - Evaluates the DAE residual written in the implicit form F(t,U,Udot)=0
764 
765   Collective
766 
767   Input Parameters:
768 + ts   - the `TS` context
769 . t    - current time
770 . U    - state vector
771 . Udot - time derivative of state vector
772 - imex - flag indicates if the method is `TSARKIMEX` so that the RHSFunction should be kept separate
773 
774   Output Parameter:
775 . Y - right-hand side
776 
777   Level: developer
778 
779   Note:
780   Most users should not need to explicitly call this routine, as it
781   is used internally within the nonlinear solvers.
782 
783   If the user did not write their equations in implicit form, this
784   function recasts them in implicit form.
785 
786 .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `TSComputeRHSFunction()`
787 @*/
788 PetscErrorCode TSComputeIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec Y, PetscBool imex)
789 {
790   TSIFunctionFn   *ifunction;
791   TSRHSFunctionFn *rhsfunction;
792   void            *ctx;
793   DM               dm;
794 
795   PetscFunctionBegin;
796   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
797   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
798   PetscValidHeaderSpecific(Udot, VEC_CLASSID, 4);
799   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
800 
801   PetscCall(TSGetDM(ts, &dm));
802   PetscCall(DMTSGetIFunction(dm, &ifunction, &ctx));
803   PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL));
804 
805   PetscCheck(rhsfunction || ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSFunction() and / or TSSetIFunction()");
806 
807   PetscCall(PetscLogEventBegin(TS_FunctionEval, U, ts, Udot, Y));
808   if (ifunction) {
809     PetscCallBack("TS callback implicit function", (*ifunction)(ts, t, U, Udot, Y, ctx));
810     ts->ifuncs++;
811   }
812   if (imex) {
813     if (!ifunction) PetscCall(VecCopy(Udot, Y));
814   } else if (rhsfunction) {
815     if (ifunction) {
816       Vec Frhs;
817 
818       PetscCall(DMGetGlobalVector(dm, &Frhs));
819       PetscCall(TSComputeRHSFunction(ts, t, U, Frhs));
820       PetscCall(VecAXPY(Y, -1, Frhs));
821       PetscCall(DMRestoreGlobalVector(dm, &Frhs));
822     } else {
823       PetscCall(TSComputeRHSFunction(ts, t, U, Y));
824       PetscCall(VecAYPX(Y, -1, Udot));
825     }
826   }
827   PetscCall(PetscLogEventEnd(TS_FunctionEval, U, ts, Udot, Y));
828   PetscFunctionReturn(PETSC_SUCCESS);
829 }
830 
831 /*
832    TSRecoverRHSJacobian - Recover the Jacobian matrix so that one can call `TSComputeRHSJacobian()` on it.
833 
834    Note:
835    This routine is needed when one switches from `TSComputeIJacobian()` to `TSComputeRHSJacobian()` because the Jacobian matrix may be shifted or scaled in `TSComputeIJacobian()`.
836 
837 */
838 static PetscErrorCode TSRecoverRHSJacobian(TS ts, Mat A, Mat B)
839 {
840   PetscFunctionBegin;
841   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
842   PetscCheck(A == ts->Arhs, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Invalid Amat");
843   PetscCheck(B == ts->Brhs, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Invalid Bmat");
844 
845   if (ts->rhsjacobian.shift) PetscCall(MatShift(A, -ts->rhsjacobian.shift));
846   if (ts->rhsjacobian.scale == -1.) PetscCall(MatScale(A, -1));
847   if (B && B == ts->Brhs && A != B) {
848     if (ts->rhsjacobian.shift) PetscCall(MatShift(B, -ts->rhsjacobian.shift));
849     if (ts->rhsjacobian.scale == -1.) PetscCall(MatScale(B, -1));
850   }
851   ts->rhsjacobian.shift = 0;
852   ts->rhsjacobian.scale = 1.;
853   PetscFunctionReturn(PETSC_SUCCESS);
854 }
855 
856 /*@
857   TSComputeIJacobian - Evaluates the Jacobian of the DAE
858 
859   Collective
860 
861   Input Parameters:
862 + ts    - the `TS` context
863 . t     - current timestep
864 . U     - state vector
865 . Udot  - time derivative of state vector
866 . shift - shift to apply, see note below
867 - imex  - flag indicates if the method is `TSARKIMEX` so that the RHSJacobian should be kept separate
868 
869   Output Parameters:
870 + A - Jacobian matrix
871 - B - matrix from which the preconditioner is constructed; often the same as `A`
872 
873   Level: developer
874 
875   Notes:
876   If F(t,U,Udot)=0 is the DAE, the required Jacobian is
877 .vb
878    dF/dU + shift*dF/dUdot
879 .ve
880   Most users should not need to explicitly call this routine, as it
881   is used internally within the nonlinear solvers.
882 
883 .seealso: [](ch_ts), `TS`, `TSSetIJacobian()`
884 @*/
885 PetscErrorCode TSComputeIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat A, Mat B, PetscBool imex)
886 {
887   TSIJacobianFn   *ijacobian;
888   TSRHSJacobianFn *rhsjacobian;
889   DM               dm;
890   void            *ctx;
891 
892   PetscFunctionBegin;
893   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
894   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
895   PetscValidHeaderSpecific(Udot, VEC_CLASSID, 4);
896   PetscValidHeaderSpecific(A, MAT_CLASSID, 6);
897   PetscValidHeaderSpecific(B, MAT_CLASSID, 7);
898 
899   PetscCall(TSGetDM(ts, &dm));
900   PetscCall(DMTSGetIJacobian(dm, &ijacobian, &ctx));
901   PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL));
902 
903   PetscCheck(rhsjacobian || ijacobian, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
904 
905   PetscCall(PetscLogEventBegin(TS_JacobianEval, U, ts, A, B));
906   if (ijacobian) {
907     PetscCallBack("TS callback implicit Jacobian", (*ijacobian)(ts, t, U, Udot, shift, A, B, ctx));
908     ts->ijacs++;
909   }
910   if (imex) {
911     if (!ijacobian) { /* system was written as Udot = G(t,U) */
912       PetscBool assembled;
913       if (rhsjacobian) {
914         Mat Arhs = NULL;
915         PetscCall(TSGetRHSMats_Private(ts, &Arhs, NULL));
916         if (A == Arhs) {
917           PetscCheck(rhsjacobian != TSComputeRHSJacobianConstant, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Unsupported operation! cannot use TSComputeRHSJacobianConstant"); /* there is no way to reconstruct shift*M-J since J cannot be reevaluated */
918           ts->rhsjacobian.time = PETSC_MIN_REAL;
919         }
920       }
921       PetscCall(MatZeroEntries(A));
922       PetscCall(MatAssembled(A, &assembled));
923       if (!assembled) {
924         PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
925         PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
926       }
927       PetscCall(MatShift(A, shift));
928       if (A != B) {
929         PetscCall(MatZeroEntries(B));
930         PetscCall(MatAssembled(B, &assembled));
931         if (!assembled) {
932           PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
933           PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
934         }
935         PetscCall(MatShift(B, shift));
936       }
937     }
938   } else {
939     Mat Arhs = NULL, Brhs = NULL;
940 
941     /* RHSJacobian needs to be converted to part of IJacobian if exists */
942     if (rhsjacobian) PetscCall(TSGetRHSMats_Private(ts, &Arhs, &Brhs));
943     if (Arhs == A) { /* No IJacobian matrix, so we only have the RHS matrix */
944       PetscObjectState Ustate;
945       PetscObjectId    Uid;
946       TSRHSFunctionFn *rhsfunction;
947 
948       PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL));
949       PetscCall(PetscObjectStateGet((PetscObject)U, &Ustate));
950       PetscCall(PetscObjectGetId((PetscObject)U, &Uid));
951       if ((rhsjacobian == TSComputeRHSJacobianConstant || (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.Xid == Uid && ts->rhsjacobian.Xstate == Ustate)) && rhsfunction != TSComputeRHSFunctionLinear)) &&
952           ts->rhsjacobian.scale == -1.) {                      /* No need to recompute RHSJacobian */
953         PetscCall(MatShift(A, shift - ts->rhsjacobian.shift)); /* revert the old shift and add the new shift with a single call to MatShift */
954         if (A != B) PetscCall(MatShift(B, shift - ts->rhsjacobian.shift));
955       } else {
956         PetscBool flg;
957 
958         if (ts->rhsjacobian.reuse) { /* Undo the damage */
959           /* MatScale has a short path for this case.
960              However, this code path is taken the first time TSComputeRHSJacobian is called
961              and the matrices have not been assembled yet */
962           PetscCall(TSRecoverRHSJacobian(ts, A, B));
963         }
964         PetscCall(TSComputeRHSJacobian(ts, t, U, A, B));
965         PetscCall(SNESGetUseMatrixFree(ts->snes, NULL, &flg));
966         /* since -snes_mf_operator uses the full SNES function it does not need to be shifted or scaled here */
967         if (!flg) {
968           PetscCall(MatScale(A, -1));
969           PetscCall(MatShift(A, shift));
970         }
971         if (A != B) {
972           PetscCall(MatScale(B, -1));
973           PetscCall(MatShift(B, shift));
974         }
975       }
976       ts->rhsjacobian.scale = -1;
977       ts->rhsjacobian.shift = shift;
978     } else if (Arhs) {  /* Both IJacobian and RHSJacobian */
979       if (!ijacobian) { /* No IJacobian provided, but we have a separate RHS matrix */
980         PetscCall(MatZeroEntries(A));
981         PetscCall(MatShift(A, shift));
982         if (A != B) {
983           PetscCall(MatZeroEntries(B));
984           PetscCall(MatShift(B, shift));
985         }
986       }
987       PetscCall(TSComputeRHSJacobian(ts, t, U, Arhs, Brhs));
988       PetscCall(MatAXPY(A, -1, Arhs, ts->axpy_pattern));
989       if (A != B) PetscCall(MatAXPY(B, -1, Brhs, ts->axpy_pattern));
990     }
991   }
992   PetscCall(PetscLogEventEnd(TS_JacobianEval, U, ts, A, B));
993   PetscFunctionReturn(PETSC_SUCCESS);
994 }
995 
996 /*@C
997   TSSetRHSFunction - Sets the routine for evaluating the function,
998   where U_t = G(t,u).
999 
1000   Logically Collective
1001 
1002   Input Parameters:
1003 + ts  - the `TS` context obtained from `TSCreate()`
1004 . r   - vector to put the computed right-hand side (or `NULL` to have it created)
1005 . f   - routine for evaluating the right-hand-side function
1006 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
1007 
1008   Level: beginner
1009 
1010   Note:
1011   You must call this function or `TSSetIFunction()` to define your ODE. You cannot use this function when solving a DAE.
1012 
1013 .seealso: [](ch_ts), `TS`, `TSRHSFunctionFn`, `TSSetRHSJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
1014 @*/
1015 PetscErrorCode TSSetRHSFunction(TS ts, Vec r, TSRHSFunctionFn *f, void *ctx)
1016 {
1017   SNES snes;
1018   Vec  ralloc = NULL;
1019   DM   dm;
1020 
1021   PetscFunctionBegin;
1022   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1023   if (r) PetscValidHeaderSpecific(r, VEC_CLASSID, 2);
1024 
1025   PetscCall(TSGetDM(ts, &dm));
1026   PetscCall(DMTSSetRHSFunction(dm, f, ctx));
1027   PetscCall(TSGetSNES(ts, &snes));
1028   if (!r && !ts->dm && ts->vec_sol) {
1029     PetscCall(VecDuplicate(ts->vec_sol, &ralloc));
1030     r = ralloc;
1031   }
1032   PetscCall(SNESSetFunction(snes, r, SNESTSFormFunction, ts));
1033   PetscCall(VecDestroy(&ralloc));
1034   PetscFunctionReturn(PETSC_SUCCESS);
1035 }
1036 
1037 /*@C
1038   TSSetSolutionFunction - Provide a function that computes the solution of the ODE or DAE
1039 
1040   Logically Collective
1041 
1042   Input Parameters:
1043 + ts  - the `TS` context obtained from `TSCreate()`
1044 . f   - routine for evaluating the solution
1045 - ctx - [optional] user-defined context for private data for the
1046           function evaluation routine (may be `NULL`)
1047 
1048   Options Database Keys:
1049 + -ts_monitor_lg_error   - create a graphical monitor of error history, requires user to have provided `TSSetSolutionFunction()`
1050 - -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
1051 
1052   Level: intermediate
1053 
1054   Notes:
1055   This routine is used for testing accuracy of time integration schemes when you already know the solution.
1056   If analytic solutions are not known for your system, consider using the Method of Manufactured Solutions to
1057   create closed-form solutions with non-physical forcing terms.
1058 
1059   For low-dimensional problems solved in serial, such as small discrete systems, `TSMonitorLGError()` can be used to monitor the error history.
1060 
1061 .seealso: [](ch_ts), `TS`, `TSSolutionFn`, `TSSetRHSJacobian()`, `TSSetIJacobian()`, `TSComputeSolutionFunction()`, `TSSetForcingFunction()`, `TSSetSolution()`, `TSGetSolution()`, `TSMonitorLGError()`, `TSMonitorDrawError()`
1062 @*/
1063 PetscErrorCode TSSetSolutionFunction(TS ts, TSSolutionFn *f, void *ctx)
1064 {
1065   DM dm;
1066 
1067   PetscFunctionBegin;
1068   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1069   PetscCall(TSGetDM(ts, &dm));
1070   PetscCall(DMTSSetSolutionFunction(dm, f, ctx));
1071   PetscFunctionReturn(PETSC_SUCCESS);
1072 }
1073 
1074 /*@C
1075   TSSetForcingFunction - Provide a function that computes a forcing term for a ODE or PDE
1076 
1077   Logically Collective
1078 
1079   Input Parameters:
1080 + ts   - the `TS` context obtained from `TSCreate()`
1081 . func - routine for evaluating the forcing function
1082 - ctx  - [optional] user-defined context for private data for the function evaluation routine
1083          (may be `NULL`)
1084 
1085   Level: intermediate
1086 
1087   Notes:
1088   This routine is useful for testing accuracy of time integration schemes when using the Method of Manufactured Solutions to
1089   create closed-form solutions with a non-physical forcing term. It allows you to use the Method of Manufactored Solution without directly editing the
1090   definition of the problem you are solving and hence possibly introducing bugs.
1091 
1092   This replaces the ODE F(u,u_t,t) = 0 the `TS` is solving with F(u,u_t,t) - func(t) = 0
1093 
1094   This forcing function does not depend on the solution to the equations, it can only depend on spatial location, time, and possibly parameters, the
1095   parameters can be passed in the ctx variable.
1096 
1097   For low-dimensional problems solved in serial, such as small discrete systems, `TSMonitorLGError()` can be used to monitor the error history.
1098 
1099 .seealso: [](ch_ts), `TS`, `TSForcingFn`, `TSSetRHSJacobian()`, `TSSetIJacobian()`,
1100 `TSComputeSolutionFunction()`, `TSSetSolutionFunction()`
1101 @*/
1102 PetscErrorCode TSSetForcingFunction(TS ts, TSForcingFn *func, void *ctx)
1103 {
1104   DM dm;
1105 
1106   PetscFunctionBegin;
1107   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1108   PetscCall(TSGetDM(ts, &dm));
1109   PetscCall(DMTSSetForcingFunction(dm, func, ctx));
1110   PetscFunctionReturn(PETSC_SUCCESS);
1111 }
1112 
1113 /*@C
1114   TSSetRHSJacobian - Sets the function to compute the Jacobian of G,
1115   where U_t = G(U,t), as well as the location to store the matrix.
1116 
1117   Logically Collective
1118 
1119   Input Parameters:
1120 + ts   - the `TS` context obtained from `TSCreate()`
1121 . Amat - (approximate) location to store Jacobian matrix entries computed by `f`
1122 . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
1123 . f    - the Jacobian evaluation routine
1124 - ctx  - [optional] user-defined context for private data for the Jacobian evaluation routine (may be `NULL`)
1125 
1126   Level: beginner
1127 
1128   Notes:
1129   You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value
1130 
1131   The `TS` solver may modify the nonzero structure and the entries of the matrices `Amat` and `Pmat` between the calls to `f()`
1132   You should not assume the values are the same in the next call to f() as you set them in the previous call.
1133 
1134 .seealso: [](ch_ts), `TS`, `TSRHSJacobianFn`, `SNESComputeJacobianDefaultColor()`,
1135 `TSSetRHSFunction()`, `TSRHSJacobianSetReuse()`, `TSSetIJacobian()`, `TSRHSFunctionFn`, `TSIFunctionFn`
1136 @*/
1137 PetscErrorCode TSSetRHSJacobian(TS ts, Mat Amat, Mat Pmat, TSRHSJacobianFn *f, void *ctx)
1138 {
1139   SNES           snes;
1140   DM             dm;
1141   TSIJacobianFn *ijacobian;
1142 
1143   PetscFunctionBegin;
1144   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1145   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
1146   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
1147   if (Amat) PetscCheckSameComm(ts, 1, Amat, 2);
1148   if (Pmat) PetscCheckSameComm(ts, 1, Pmat, 3);
1149 
1150   PetscCall(TSGetDM(ts, &dm));
1151   PetscCall(DMTSSetRHSJacobian(dm, f, ctx));
1152   PetscCall(DMTSGetIJacobian(dm, &ijacobian, NULL));
1153   PetscCall(TSGetSNES(ts, &snes));
1154   if (!ijacobian) PetscCall(SNESSetJacobian(snes, Amat, Pmat, SNESTSFormJacobian, ts));
1155   if (Amat) {
1156     PetscCall(PetscObjectReference((PetscObject)Amat));
1157     PetscCall(MatDestroy(&ts->Arhs));
1158     ts->Arhs = Amat;
1159   }
1160   if (Pmat) {
1161     PetscCall(PetscObjectReference((PetscObject)Pmat));
1162     PetscCall(MatDestroy(&ts->Brhs));
1163     ts->Brhs = Pmat;
1164   }
1165   PetscFunctionReturn(PETSC_SUCCESS);
1166 }
1167 
1168 /*@C
1169   TSSetIFunction - Set the function to compute F(t,U,U_t) where F() = 0 is the DAE to be solved.
1170 
1171   Logically Collective
1172 
1173   Input Parameters:
1174 + ts  - the `TS` context obtained from `TSCreate()`
1175 . r   - vector to hold the residual (or `NULL` to have it created internally)
1176 . f   - the function evaluation routine
1177 - ctx - user-defined context for private data for the function evaluation routine (may be `NULL`)
1178 
1179   Level: beginner
1180 
1181   Note:
1182   The user MUST call either this routine or `TSSetRHSFunction()` to define the ODE.  When solving DAEs you must use this function.
1183 
1184 .seealso: [](ch_ts), `TS`, `TSIFunctionFn`, `TSSetRHSJacobian()`, `TSSetRHSFunction()`,
1185 `TSSetIJacobian()`
1186 @*/
1187 PetscErrorCode TSSetIFunction(TS ts, Vec r, TSIFunctionFn *f, void *ctx)
1188 {
1189   SNES snes;
1190   Vec  ralloc = NULL;
1191   DM   dm;
1192 
1193   PetscFunctionBegin;
1194   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1195   if (r) PetscValidHeaderSpecific(r, VEC_CLASSID, 2);
1196 
1197   PetscCall(TSGetDM(ts, &dm));
1198   PetscCall(DMTSSetIFunction(dm, f, ctx));
1199 
1200   PetscCall(TSGetSNES(ts, &snes));
1201   if (!r && !ts->dm && ts->vec_sol) {
1202     PetscCall(VecDuplicate(ts->vec_sol, &ralloc));
1203     r = ralloc;
1204   }
1205   PetscCall(SNESSetFunction(snes, r, SNESTSFormFunction, ts));
1206   PetscCall(VecDestroy(&ralloc));
1207   PetscFunctionReturn(PETSC_SUCCESS);
1208 }
1209 
1210 /*@C
1211   TSGetIFunction - Returns the vector where the implicit residual is stored and the function/context to compute it.
1212 
1213   Not Collective
1214 
1215   Input Parameter:
1216 . ts - the `TS` context
1217 
1218   Output Parameters:
1219 + r    - vector to hold residual (or `NULL`)
1220 . func - the function to compute residual (or `NULL`)
1221 - ctx  - the function context (or `NULL`)
1222 
1223   Level: advanced
1224 
1225 .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `SNESGetFunction()`
1226 @*/
1227 PetscErrorCode TSGetIFunction(TS ts, Vec *r, TSIFunctionFn **func, void **ctx)
1228 {
1229   SNES snes;
1230   DM   dm;
1231 
1232   PetscFunctionBegin;
1233   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1234   PetscCall(TSGetSNES(ts, &snes));
1235   PetscCall(SNESGetFunction(snes, r, NULL, NULL));
1236   PetscCall(TSGetDM(ts, &dm));
1237   PetscCall(DMTSGetIFunction(dm, func, ctx));
1238   PetscFunctionReturn(PETSC_SUCCESS);
1239 }
1240 
1241 /*@C
1242   TSGetRHSFunction - Returns the vector where the right-hand side is stored and the function/context to compute it.
1243 
1244   Not Collective
1245 
1246   Input Parameter:
1247 . ts - the `TS` context
1248 
1249   Output Parameters:
1250 + r    - vector to hold computed right-hand side (or `NULL`)
1251 . func - the function to compute right-hand side (or `NULL`)
1252 - ctx  - the function context (or `NULL`)
1253 
1254   Level: advanced
1255 
1256 .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `SNESGetFunction()`
1257 @*/
1258 PetscErrorCode TSGetRHSFunction(TS ts, Vec *r, TSRHSFunctionFn **func, void **ctx)
1259 {
1260   SNES snes;
1261   DM   dm;
1262 
1263   PetscFunctionBegin;
1264   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1265   PetscCall(TSGetSNES(ts, &snes));
1266   PetscCall(SNESGetFunction(snes, r, NULL, NULL));
1267   PetscCall(TSGetDM(ts, &dm));
1268   PetscCall(DMTSGetRHSFunction(dm, func, ctx));
1269   PetscFunctionReturn(PETSC_SUCCESS);
1270 }
1271 
1272 /*@C
1273   TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
1274   provided with `TSSetIFunction()`.
1275 
1276   Logically Collective
1277 
1278   Input Parameters:
1279 + ts   - the `TS` context obtained from `TSCreate()`
1280 . Amat - (approximate) matrix to store Jacobian entries computed by `f`
1281 . Pmat - matrix used to compute preconditioner (usually the same as `Amat`)
1282 . f    - the Jacobian evaluation routine
1283 - ctx  - user-defined context for private data for the Jacobian evaluation routine (may be `NULL`)
1284 
1285   Level: beginner
1286 
1287   Notes:
1288   The matrices `Amat` and `Pmat` are exactly the matrices that are used by `SNES` for the nonlinear solve.
1289 
1290   If you know the operator Amat has a null space you can use `MatSetNullSpace()` and `MatSetTransposeNullSpace()` to supply the null
1291   space to `Amat` and the `KSP` solvers will automatically use that null space as needed during the solution process.
1292 
1293   The matrix dF/dU + a*dF/dU_t you provide turns out to be
1294   the Jacobian of F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
1295   The time integrator internally approximates U_t by W+a*U where the positive "shift"
1296   a and vector W depend on the integration method, step size, and past states. For example with
1297   the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
1298   W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
1299 
1300   You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value
1301 
1302   The TS solver may modify the nonzero structure and the entries of the matrices `Amat` and `Pmat` between the calls to `f`
1303   You should not assume the values are the same in the next call to `f` as you set them in the previous call.
1304 
1305   In case `TSSetRHSJacobian()` is also used in conjunction with a fully-implicit solver,
1306   multilevel linear solvers, e.g. `PCMG`, will likely not work due to the way `TS` handles rhs matrices.
1307 
1308 .seealso: [](ch_ts), `TS`, `TSIJacobianFn`, `TSSetIFunction()`, `TSSetRHSJacobian()`,
1309 `SNESComputeJacobianDefaultColor()`, `SNESComputeJacobianDefault()`, `TSSetRHSFunction()`
1310 @*/
1311 PetscErrorCode TSSetIJacobian(TS ts, Mat Amat, Mat Pmat, TSIJacobianFn *f, void *ctx)
1312 {
1313   SNES snes;
1314   DM   dm;
1315 
1316   PetscFunctionBegin;
1317   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1318   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
1319   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
1320   if (Amat) PetscCheckSameComm(ts, 1, Amat, 2);
1321   if (Pmat) PetscCheckSameComm(ts, 1, Pmat, 3);
1322 
1323   PetscCall(TSGetDM(ts, &dm));
1324   PetscCall(DMTSSetIJacobian(dm, f, ctx));
1325 
1326   PetscCall(TSGetSNES(ts, &snes));
1327   PetscCall(SNESSetJacobian(snes, Amat, Pmat, SNESTSFormJacobian, ts));
1328   PetscFunctionReturn(PETSC_SUCCESS);
1329 }
1330 
1331 /*@
1332   TSRHSJacobianSetReuse - restore the RHS Jacobian before calling the user-provided `TSRHSJacobianFn` function again
1333 
1334   Logically Collective
1335 
1336   Input Parameters:
1337 + ts    - `TS` context obtained from `TSCreate()`
1338 - reuse - `PETSC_TRUE` if the RHS Jacobian
1339 
1340   Level: intermediate
1341 
1342   Notes:
1343   Without this flag, `TS` will change the sign and shift the RHS Jacobian for a
1344   finite-time-step implicit solve, in which case the user function will need to recompute the
1345   entire Jacobian.  The `reuse `flag must be set if the evaluation function assumes that the
1346   matrix entries have not been changed by the `TS`.
1347 
1348 .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `TSComputeRHSJacobianConstant()`
1349 @*/
1350 PetscErrorCode TSRHSJacobianSetReuse(TS ts, PetscBool reuse)
1351 {
1352   PetscFunctionBegin;
1353   ts->rhsjacobian.reuse = reuse;
1354   PetscFunctionReturn(PETSC_SUCCESS);
1355 }
1356 
1357 /*@C
1358   TSSetI2Function - Set the function to compute F(t,U,U_t,U_tt) where F = 0 is the DAE to be solved.
1359 
1360   Logically Collective
1361 
1362   Input Parameters:
1363 + ts  - the `TS` context obtained from `TSCreate()`
1364 . F   - vector to hold the residual (or `NULL` to have it created internally)
1365 . fun - the function evaluation routine
1366 - ctx - user-defined context for private data for the function evaluation routine (may be `NULL`)
1367 
1368   Level: beginner
1369 
1370 .seealso: [](ch_ts), `TS`, `TSI2FunctionFn`, `TSSetI2Jacobian()`, `TSSetIFunction()`,
1371 `TSCreate()`, `TSSetRHSFunction()`
1372 @*/
1373 PetscErrorCode TSSetI2Function(TS ts, Vec F, TSI2FunctionFn *fun, void *ctx)
1374 {
1375   DM dm;
1376 
1377   PetscFunctionBegin;
1378   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1379   if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 2);
1380   PetscCall(TSSetIFunction(ts, F, NULL, NULL));
1381   PetscCall(TSGetDM(ts, &dm));
1382   PetscCall(DMTSSetI2Function(dm, fun, ctx));
1383   PetscFunctionReturn(PETSC_SUCCESS);
1384 }
1385 
1386 /*@C
1387   TSGetI2Function - Returns the vector where the implicit residual is stored and the function/context to compute it.
1388 
1389   Not Collective
1390 
1391   Input Parameter:
1392 . ts - the `TS` context
1393 
1394   Output Parameters:
1395 + r   - vector to hold residual (or `NULL`)
1396 . fun - the function to compute residual (or `NULL`)
1397 - ctx - the function context (or `NULL`)
1398 
1399   Level: advanced
1400 
1401 .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `SNESGetFunction()`, `TSCreate()`
1402 @*/
1403 PetscErrorCode TSGetI2Function(TS ts, Vec *r, TSI2FunctionFn **fun, void **ctx)
1404 {
1405   SNES snes;
1406   DM   dm;
1407 
1408   PetscFunctionBegin;
1409   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1410   PetscCall(TSGetSNES(ts, &snes));
1411   PetscCall(SNESGetFunction(snes, r, NULL, NULL));
1412   PetscCall(TSGetDM(ts, &dm));
1413   PetscCall(DMTSGetI2Function(dm, fun, ctx));
1414   PetscFunctionReturn(PETSC_SUCCESS);
1415 }
1416 
1417 /*@C
1418   TSSetI2Jacobian - Set the function to compute the matrix dF/dU + v*dF/dU_t  + a*dF/dU_tt
1419   where F(t,U,U_t,U_tt) is the function you provided with `TSSetI2Function()`.
1420 
1421   Logically Collective
1422 
1423   Input Parameters:
1424 + ts  - the `TS` context obtained from `TSCreate()`
1425 . J   - matrix to hold the Jacobian values
1426 . P   - matrix for constructing the preconditioner (may be same as `J`)
1427 . jac - the Jacobian evaluation routine, see `TSI2JacobianFn` for the calling sequence
1428 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be `NULL`)
1429 
1430   Level: beginner
1431 
1432   Notes:
1433   The matrices `J` and `P` are exactly the matrices that are used by `SNES` for the nonlinear solve.
1434 
1435   The matrix dF/dU + v*dF/dU_t + a*dF/dU_tt you provide turns out to be
1436   the Jacobian of G(U) = F(t,U,W+v*U,W'+a*U) where F(t,U,U_t,U_tt) = 0 is the DAE to be solved.
1437   The time integrator internally approximates U_t by W+v*U and U_tt by W'+a*U  where the positive "shift"
1438   parameters 'v' and 'a' and vectors W, W' depend on the integration method, step size, and past states.
1439 
1440 .seealso: [](ch_ts), `TS`, `TSI2JacobianFn`, `TSSetI2Function()`, `TSGetI2Jacobian()`
1441 @*/
1442 PetscErrorCode TSSetI2Jacobian(TS ts, Mat J, Mat P, TSI2JacobianFn *jac, void *ctx)
1443 {
1444   DM dm;
1445 
1446   PetscFunctionBegin;
1447   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1448   if (J) PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
1449   if (P) PetscValidHeaderSpecific(P, MAT_CLASSID, 3);
1450   PetscCall(TSSetIJacobian(ts, J, P, NULL, NULL));
1451   PetscCall(TSGetDM(ts, &dm));
1452   PetscCall(DMTSSetI2Jacobian(dm, jac, ctx));
1453   PetscFunctionReturn(PETSC_SUCCESS);
1454 }
1455 
1456 /*@C
1457   TSGetI2Jacobian - Returns the implicit Jacobian at the present timestep.
1458 
1459   Not Collective, but parallel objects are returned if `TS` is parallel
1460 
1461   Input Parameter:
1462 . ts - The `TS` context obtained from `TSCreate()`
1463 
1464   Output Parameters:
1465 + J   - The (approximate) Jacobian of F(t,U,U_t,U_tt)
1466 . P   - The matrix from which the preconditioner is constructed, often the same as `J`
1467 . jac - The function to compute the Jacobian matrices
1468 - ctx - User-defined context for Jacobian evaluation routine
1469 
1470   Level: advanced
1471 
1472   Note:
1473   You can pass in `NULL` for any return argument you do not need.
1474 
1475 .seealso: [](ch_ts), `TS`, `TSGetTimeStep()`, `TSGetMatrices()`, `TSGetTime()`, `TSGetStepNumber()`, `TSSetI2Jacobian()`, `TSGetI2Function()`, `TSCreate()`
1476 @*/
1477 PetscErrorCode TSGetI2Jacobian(TS ts, Mat *J, Mat *P, TSI2JacobianFn **jac, void **ctx)
1478 {
1479   SNES snes;
1480   DM   dm;
1481 
1482   PetscFunctionBegin;
1483   PetscCall(TSGetSNES(ts, &snes));
1484   PetscCall(SNESSetUpMatrices(snes));
1485   PetscCall(SNESGetJacobian(snes, J, P, NULL, NULL));
1486   PetscCall(TSGetDM(ts, &dm));
1487   PetscCall(DMTSGetI2Jacobian(dm, jac, ctx));
1488   PetscFunctionReturn(PETSC_SUCCESS);
1489 }
1490 
1491 /*@
1492   TSComputeI2Function - Evaluates the DAE residual written in implicit form F(t,U,U_t,U_tt) = 0
1493 
1494   Collective
1495 
1496   Input Parameters:
1497 + ts - the `TS` context
1498 . t  - current time
1499 . U  - state vector
1500 . V  - time derivative of state vector (U_t)
1501 - A  - second time derivative of state vector (U_tt)
1502 
1503   Output Parameter:
1504 . F - the residual vector
1505 
1506   Level: developer
1507 
1508   Note:
1509   Most users should not need to explicitly call this routine, as it
1510   is used internally within the nonlinear solvers.
1511 
1512 .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `TSGetI2Function()`
1513 @*/
1514 PetscErrorCode TSComputeI2Function(TS ts, PetscReal t, Vec U, Vec V, Vec A, Vec F)
1515 {
1516   DM               dm;
1517   TSI2FunctionFn  *I2Function;
1518   void            *ctx;
1519   TSRHSFunctionFn *rhsfunction;
1520 
1521   PetscFunctionBegin;
1522   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1523   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
1524   PetscValidHeaderSpecific(V, VEC_CLASSID, 4);
1525   PetscValidHeaderSpecific(A, VEC_CLASSID, 5);
1526   PetscValidHeaderSpecific(F, VEC_CLASSID, 6);
1527 
1528   PetscCall(TSGetDM(ts, &dm));
1529   PetscCall(DMTSGetI2Function(dm, &I2Function, &ctx));
1530   PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL));
1531 
1532   if (!I2Function) {
1533     PetscCall(TSComputeIFunction(ts, t, U, A, F, PETSC_FALSE));
1534     PetscFunctionReturn(PETSC_SUCCESS);
1535   }
1536 
1537   PetscCall(PetscLogEventBegin(TS_FunctionEval, U, ts, V, F));
1538 
1539   PetscCallBack("TS callback implicit function", I2Function(ts, t, U, V, A, F, ctx));
1540 
1541   if (rhsfunction) {
1542     Vec Frhs;
1543 
1544     PetscCall(DMGetGlobalVector(dm, &Frhs));
1545     PetscCall(TSComputeRHSFunction(ts, t, U, Frhs));
1546     PetscCall(VecAXPY(F, -1, Frhs));
1547     PetscCall(DMRestoreGlobalVector(dm, &Frhs));
1548   }
1549 
1550   PetscCall(PetscLogEventEnd(TS_FunctionEval, U, ts, V, F));
1551   PetscFunctionReturn(PETSC_SUCCESS);
1552 }
1553 
1554 /*@
1555   TSComputeI2Jacobian - Evaluates the Jacobian of the DAE
1556 
1557   Collective
1558 
1559   Input Parameters:
1560 + ts     - the `TS` context
1561 . t      - current timestep
1562 . U      - state vector
1563 . V      - time derivative of state vector
1564 . A      - second time derivative of state vector
1565 . shiftV - shift to apply, see note below
1566 - shiftA - shift to apply, see note below
1567 
1568   Output Parameters:
1569 + J - Jacobian matrix
1570 - P - optional preconditioning matrix
1571 
1572   Level: developer
1573 
1574   Notes:
1575   If F(t,U,V,A)=0 is the DAE, the required Jacobian is
1576 
1577   dF/dU + shiftV*dF/dV + shiftA*dF/dA
1578 
1579   Most users should not need to explicitly call this routine, as it
1580   is used internally within the nonlinear solvers.
1581 
1582 .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`
1583 @*/
1584 PetscErrorCode TSComputeI2Jacobian(TS ts, PetscReal t, Vec U, Vec V, Vec A, PetscReal shiftV, PetscReal shiftA, Mat J, Mat P)
1585 {
1586   DM               dm;
1587   TSI2JacobianFn  *I2Jacobian;
1588   void            *ctx;
1589   TSRHSJacobianFn *rhsjacobian;
1590 
1591   PetscFunctionBegin;
1592   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1593   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
1594   PetscValidHeaderSpecific(V, VEC_CLASSID, 4);
1595   PetscValidHeaderSpecific(A, VEC_CLASSID, 5);
1596   PetscValidHeaderSpecific(J, MAT_CLASSID, 8);
1597   PetscValidHeaderSpecific(P, MAT_CLASSID, 9);
1598 
1599   PetscCall(TSGetDM(ts, &dm));
1600   PetscCall(DMTSGetI2Jacobian(dm, &I2Jacobian, &ctx));
1601   PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL));
1602 
1603   if (!I2Jacobian) {
1604     PetscCall(TSComputeIJacobian(ts, t, U, A, shiftA, J, P, PETSC_FALSE));
1605     PetscFunctionReturn(PETSC_SUCCESS);
1606   }
1607 
1608   PetscCall(PetscLogEventBegin(TS_JacobianEval, U, ts, J, P));
1609   PetscCallBack("TS callback implicit Jacobian", I2Jacobian(ts, t, U, V, A, shiftV, shiftA, J, P, ctx));
1610   if (rhsjacobian) {
1611     Mat Jrhs, Prhs;
1612     PetscCall(TSGetRHSMats_Private(ts, &Jrhs, &Prhs));
1613     PetscCall(TSComputeRHSJacobian(ts, t, U, Jrhs, Prhs));
1614     PetscCall(MatAXPY(J, -1, Jrhs, ts->axpy_pattern));
1615     if (P != J) PetscCall(MatAXPY(P, -1, Prhs, ts->axpy_pattern));
1616   }
1617 
1618   PetscCall(PetscLogEventEnd(TS_JacobianEval, U, ts, J, P));
1619   PetscFunctionReturn(PETSC_SUCCESS);
1620 }
1621 
1622 /*@C
1623   TSSetTransientVariable - sets function to transform from state to transient variables
1624 
1625   Logically Collective
1626 
1627   Input Parameters:
1628 + ts   - time stepping context on which to change the transient variable
1629 . tvar - a function that transforms to transient variables, see `TSTransientVariableFn` for the calling sequence
1630 - ctx  - a context for tvar
1631 
1632   Level: advanced
1633 
1634   Notes:
1635   This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., `TSBDF`)
1636   can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
1637   well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
1638   C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
1639   evaluated via the chain rule, as in
1640 .vb
1641      dF/dP + shift * dF/dCdot dC/dP.
1642 .ve
1643 
1644 .seealso: [](ch_ts), `TS`, `TSBDF`, `TSTransientVariableFn`, `DMTSSetTransientVariable()`, `DMTSGetTransientVariable()`, `TSSetIFunction()`, `TSSetIJacobian()`
1645 @*/
1646 PetscErrorCode TSSetTransientVariable(TS ts, TSTransientVariableFn *tvar, void *ctx)
1647 {
1648   DM dm;
1649 
1650   PetscFunctionBegin;
1651   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1652   PetscCall(TSGetDM(ts, &dm));
1653   PetscCall(DMTSSetTransientVariable(dm, tvar, ctx));
1654   PetscFunctionReturn(PETSC_SUCCESS);
1655 }
1656 
1657 /*@
1658   TSComputeTransientVariable - transforms state (primitive) variables to transient (conservative) variables
1659 
1660   Logically Collective
1661 
1662   Input Parameters:
1663 + ts - TS on which to compute
1664 - U  - state vector to be transformed to transient variables
1665 
1666   Output Parameter:
1667 . C - transient (conservative) variable
1668 
1669   Level: developer
1670 
1671   Developer Notes:
1672   If `DMTSSetTransientVariable()` has not been called, then C is not modified in this routine and C = `NULL` is allowed.
1673   This makes it safe to call without a guard.  One can use `TSHasTransientVariable()` to check if transient variables are
1674   being used.
1675 
1676 .seealso: [](ch_ts), `TS`, `TSBDF`, `DMTSSetTransientVariable()`, `TSComputeIFunction()`, `TSComputeIJacobian()`
1677 @*/
1678 PetscErrorCode TSComputeTransientVariable(TS ts, Vec U, Vec C)
1679 {
1680   DM   dm;
1681   DMTS dmts;
1682 
1683   PetscFunctionBegin;
1684   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1685   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
1686   PetscCall(TSGetDM(ts, &dm));
1687   PetscCall(DMGetDMTS(dm, &dmts));
1688   if (dmts->ops->transientvar) {
1689     PetscValidHeaderSpecific(C, VEC_CLASSID, 3);
1690     PetscCall((*dmts->ops->transientvar)(ts, U, C, dmts->transientvarctx));
1691   }
1692   PetscFunctionReturn(PETSC_SUCCESS);
1693 }
1694 
1695 /*@
1696   TSHasTransientVariable - determine whether transient variables have been set
1697 
1698   Logically Collective
1699 
1700   Input Parameter:
1701 . ts - `TS` on which to compute
1702 
1703   Output Parameter:
1704 . has - `PETSC_TRUE` if transient variables have been set
1705 
1706   Level: developer
1707 
1708 .seealso: [](ch_ts), `TS`, `TSBDF`, `DMTSSetTransientVariable()`, `TSComputeTransientVariable()`
1709 @*/
1710 PetscErrorCode TSHasTransientVariable(TS ts, PetscBool *has)
1711 {
1712   DM   dm;
1713   DMTS dmts;
1714 
1715   PetscFunctionBegin;
1716   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1717   PetscCall(TSGetDM(ts, &dm));
1718   PetscCall(DMGetDMTS(dm, &dmts));
1719   *has = dmts->ops->transientvar ? PETSC_TRUE : PETSC_FALSE;
1720   PetscFunctionReturn(PETSC_SUCCESS);
1721 }
1722 
1723 /*@
1724   TS2SetSolution - Sets the initial solution and time derivative vectors
1725   for use by the `TS` routines handling second order equations.
1726 
1727   Logically Collective
1728 
1729   Input Parameters:
1730 + ts - the `TS` context obtained from `TSCreate()`
1731 . u  - the solution vector
1732 - v  - the time derivative vector
1733 
1734   Level: beginner
1735 
1736 .seealso: [](ch_ts), `TS`
1737 @*/
1738 PetscErrorCode TS2SetSolution(TS ts, Vec u, Vec v)
1739 {
1740   PetscFunctionBegin;
1741   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1742   PetscValidHeaderSpecific(u, VEC_CLASSID, 2);
1743   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
1744   PetscCall(TSSetSolution(ts, u));
1745   PetscCall(PetscObjectReference((PetscObject)v));
1746   PetscCall(VecDestroy(&ts->vec_dot));
1747   ts->vec_dot = v;
1748   PetscFunctionReturn(PETSC_SUCCESS);
1749 }
1750 
1751 /*@
1752   TS2GetSolution - Returns the solution and time derivative at the present timestep
1753   for second order equations.
1754 
1755   Not Collective
1756 
1757   Input Parameter:
1758 . ts - the `TS` context obtained from `TSCreate()`
1759 
1760   Output Parameters:
1761 + u - the vector containing the solution
1762 - v - the vector containing the time derivative
1763 
1764   Level: intermediate
1765 
1766   Notes:
1767   It is valid to call this routine inside the function
1768   that you are evaluating in order to move to the new timestep. This vector not
1769   changed until the solution at the next timestep has been calculated.
1770 
1771 .seealso: [](ch_ts), `TS`, `TS2SetSolution()`, `TSGetTimeStep()`, `TSGetTime()`
1772 @*/
1773 PetscErrorCode TS2GetSolution(TS ts, Vec *u, Vec *v)
1774 {
1775   PetscFunctionBegin;
1776   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1777   if (u) PetscAssertPointer(u, 2);
1778   if (v) PetscAssertPointer(v, 3);
1779   if (u) *u = ts->vec_sol;
1780   if (v) *v = ts->vec_dot;
1781   PetscFunctionReturn(PETSC_SUCCESS);
1782 }
1783 
1784 /*@
1785   TSLoad - Loads a `TS` that has been stored in binary  with `TSView()`.
1786 
1787   Collective
1788 
1789   Input Parameters:
1790 + ts     - the newly loaded `TS`, this needs to have been created with `TSCreate()` or
1791            some related function before a call to `TSLoad()`.
1792 - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`
1793 
1794   Level: intermediate
1795 
1796   Note:
1797   The type is determined by the data in the file, any type set into the `TS` before this call is ignored.
1798 
1799 .seealso: [](ch_ts), `TS`, `PetscViewer`, `PetscViewerBinaryOpen()`, `TSView()`, `MatLoad()`, `VecLoad()`
1800 @*/
1801 PetscErrorCode TSLoad(TS ts, PetscViewer viewer)
1802 {
1803   PetscBool isbinary;
1804   PetscInt  classid;
1805   char      type[256];
1806   DMTS      sdm;
1807   DM        dm;
1808 
1809   PetscFunctionBegin;
1810   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1811   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1812   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1813   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1814 
1815   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
1816   PetscCheck(classid == TS_FILE_CLASSID, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Not TS next in file");
1817   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
1818   PetscCall(TSSetType(ts, type));
1819   PetscTryTypeMethod(ts, load, viewer);
1820   PetscCall(DMCreate(PetscObjectComm((PetscObject)ts), &dm));
1821   PetscCall(DMLoad(dm, viewer));
1822   PetscCall(TSSetDM(ts, dm));
1823   PetscCall(DMCreateGlobalVector(ts->dm, &ts->vec_sol));
1824   PetscCall(VecLoad(ts->vec_sol, viewer));
1825   PetscCall(DMGetDMTS(ts->dm, &sdm));
1826   PetscCall(DMTSLoad(sdm, viewer));
1827   PetscFunctionReturn(PETSC_SUCCESS);
1828 }
1829 
1830 #include <petscdraw.h>
1831 #if defined(PETSC_HAVE_SAWS)
1832   #include <petscviewersaws.h>
1833 #endif
1834 
1835 /*@
1836   TSViewFromOptions - View a `TS` based on values in the options database
1837 
1838   Collective
1839 
1840   Input Parameters:
1841 + ts   - the `TS` context
1842 . obj  - Optional object that provides the prefix for the options database keys
1843 - name - command line option string to be passed by user
1844 
1845   Level: intermediate
1846 
1847 .seealso: [](ch_ts), `TS`, `TSView`, `PetscObjectViewFromOptions()`, `TSCreate()`
1848 @*/
1849 PetscErrorCode TSViewFromOptions(TS ts, PetscObject obj, const char name[])
1850 {
1851   PetscFunctionBegin;
1852   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1853   PetscCall(PetscObjectViewFromOptions((PetscObject)ts, obj, name));
1854   PetscFunctionReturn(PETSC_SUCCESS);
1855 }
1856 
1857 /*@
1858   TSView - Prints the `TS` data structure.
1859 
1860   Collective
1861 
1862   Input Parameters:
1863 + ts     - the `TS` context obtained from `TSCreate()`
1864 - viewer - visualization context
1865 
1866   Options Database Key:
1867 . -ts_view - calls `TSView()` at end of `TSStep()`
1868 
1869   Level: beginner
1870 
1871   Notes:
1872   The available visualization contexts include
1873 +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1874 -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
1875   output where only the first processor opens
1876   the file.  All other processors send their
1877   data to the first processor to print.
1878 
1879   The user can open an alternative visualization context with
1880   `PetscViewerASCIIOpen()` - output to a specified file.
1881 
1882   In the debugger you can do call `TSView`(ts,0) to display the `TS` solver. (The same holds for any PETSc object viewer).
1883 
1884 .seealso: [](ch_ts), `TS`, `PetscViewer`, `PetscViewerASCIIOpen()`
1885 @*/
1886 PetscErrorCode TSView(TS ts, PetscViewer viewer)
1887 {
1888   TSType    type;
1889   PetscBool iascii, isstring, isundials, isbinary, isdraw;
1890   DMTS      sdm;
1891 #if defined(PETSC_HAVE_SAWS)
1892   PetscBool issaws;
1893 #endif
1894 
1895   PetscFunctionBegin;
1896   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1897   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts), &viewer));
1898   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1899   PetscCheckSameComm(ts, 1, viewer, 2);
1900 
1901   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1902   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1903   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1904   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1905 #if defined(PETSC_HAVE_SAWS)
1906   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
1907 #endif
1908   if (iascii) {
1909     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ts, viewer));
1910     if (ts->ops->view) {
1911       PetscCall(PetscViewerASCIIPushTab(viewer));
1912       PetscUseTypeMethod(ts, view, viewer);
1913       PetscCall(PetscViewerASCIIPopTab(viewer));
1914     }
1915     if (ts->max_steps < PETSC_INT_MAX) PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum steps=%" PetscInt_FMT "\n", ts->max_steps));
1916     if (ts->max_time < PETSC_MAX_REAL) PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum time=%g\n", (double)ts->max_time));
1917     if (ts->ifuncs) PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of I function evaluations=%" PetscInt_FMT "\n", ts->ifuncs));
1918     if (ts->ijacs) PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of I Jacobian evaluations=%" PetscInt_FMT "\n", ts->ijacs));
1919     if (ts->rhsfuncs) PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of RHS function evaluations=%" PetscInt_FMT "\n", ts->rhsfuncs));
1920     if (ts->rhsjacs) PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of RHS Jacobian evaluations=%" PetscInt_FMT "\n", ts->rhsjacs));
1921     if (ts->usessnes) {
1922       PetscBool lin;
1923       if (ts->problem_type == TS_NONLINEAR) PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of nonlinear solver iterations=%" PetscInt_FMT "\n", ts->snes_its));
1924       PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of linear solver iterations=%" PetscInt_FMT "\n", ts->ksp_its));
1925       PetscCall(PetscObjectTypeCompareAny((PetscObject)ts->snes, &lin, SNESKSPONLY, SNESKSPTRANSPOSEONLY, ""));
1926       PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of %slinear solve failures=%" PetscInt_FMT "\n", lin ? "" : "non", ts->num_snes_failures));
1927     }
1928     PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of rejected steps=%" PetscInt_FMT "\n", ts->reject));
1929     if (ts->vrtol) PetscCall(PetscViewerASCIIPrintf(viewer, "  using vector of relative error tolerances, "));
1930     else PetscCall(PetscViewerASCIIPrintf(viewer, "  using relative error tolerance of %g, ", (double)ts->rtol));
1931     if (ts->vatol) PetscCall(PetscViewerASCIIPrintf(viewer, "  using vector of absolute error tolerances\n"));
1932     else PetscCall(PetscViewerASCIIPrintf(viewer, "  using absolute error tolerance of %g\n", (double)ts->atol));
1933     PetscCall(PetscViewerASCIIPushTab(viewer));
1934     PetscCall(TSAdaptView(ts->adapt, viewer));
1935     PetscCall(PetscViewerASCIIPopTab(viewer));
1936   } else if (isstring) {
1937     PetscCall(TSGetType(ts, &type));
1938     PetscCall(PetscViewerStringSPrintf(viewer, " TSType: %-7.7s", type));
1939     PetscTryTypeMethod(ts, view, viewer);
1940   } else if (isbinary) {
1941     PetscInt    classid = TS_FILE_CLASSID;
1942     MPI_Comm    comm;
1943     PetscMPIInt rank;
1944     char        type[256];
1945 
1946     PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
1947     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1948     if (rank == 0) {
1949       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
1950       PetscCall(PetscStrncpy(type, ((PetscObject)ts)->type_name, 256));
1951       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
1952     }
1953     PetscTryTypeMethod(ts, view, viewer);
1954     if (ts->adapt) PetscCall(TSAdaptView(ts->adapt, viewer));
1955     PetscCall(DMView(ts->dm, viewer));
1956     PetscCall(VecView(ts->vec_sol, viewer));
1957     PetscCall(DMGetDMTS(ts->dm, &sdm));
1958     PetscCall(DMTSView(sdm, viewer));
1959   } else if (isdraw) {
1960     PetscDraw draw;
1961     char      str[36];
1962     PetscReal x, y, bottom, h;
1963 
1964     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1965     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
1966     PetscCall(PetscStrncpy(str, "TS: ", sizeof(str)));
1967     PetscCall(PetscStrlcat(str, ((PetscObject)ts)->type_name, sizeof(str)));
1968     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_BLACK, PETSC_DRAW_BLACK, str, NULL, &h));
1969     bottom = y - h;
1970     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1971     PetscTryTypeMethod(ts, view, viewer);
1972     if (ts->adapt) PetscCall(TSAdaptView(ts->adapt, viewer));
1973     if (ts->snes) PetscCall(SNESView(ts->snes, viewer));
1974     PetscCall(PetscDrawPopCurrentPoint(draw));
1975 #if defined(PETSC_HAVE_SAWS)
1976   } else if (issaws) {
1977     PetscMPIInt rank;
1978     const char *name;
1979 
1980     PetscCall(PetscObjectGetName((PetscObject)ts, &name));
1981     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1982     if (!((PetscObject)ts)->amsmem && rank == 0) {
1983       char dir[1024];
1984 
1985       PetscCall(PetscObjectViewSAWs((PetscObject)ts, viewer));
1986       PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Objects/%s/time_step", name));
1987       PetscCallSAWs(SAWs_Register, (dir, &ts->steps, 1, SAWs_READ, SAWs_INT));
1988       PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Objects/%s/time", name));
1989       PetscCallSAWs(SAWs_Register, (dir, &ts->ptime, 1, SAWs_READ, SAWs_DOUBLE));
1990     }
1991     PetscTryTypeMethod(ts, view, viewer);
1992 #endif
1993   }
1994   if (ts->snes && ts->usessnes) {
1995     PetscCall(PetscViewerASCIIPushTab(viewer));
1996     PetscCall(SNESView(ts->snes, viewer));
1997     PetscCall(PetscViewerASCIIPopTab(viewer));
1998   }
1999   PetscCall(DMGetDMTS(ts->dm, &sdm));
2000   PetscCall(DMTSView(sdm, viewer));
2001 
2002   PetscCall(PetscViewerASCIIPushTab(viewer));
2003   PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSSUNDIALS, &isundials));
2004   PetscCall(PetscViewerASCIIPopTab(viewer));
2005   PetscFunctionReturn(PETSC_SUCCESS);
2006 }
2007 
2008 /*@
2009   TSSetApplicationContext - Sets an optional user-defined context for
2010   the timesteppers.
2011 
2012   Logically Collective
2013 
2014   Input Parameters:
2015 + ts  - the `TS` context obtained from `TSCreate()`
2016 - ctx - user context
2017 
2018   Level: intermediate
2019 
2020   Fortran Notes:
2021   You must write a Fortran interface definition for this
2022   function that tells Fortran the Fortran derived data type that you are passing in as the `ctx` argument.
2023 
2024 .seealso: [](ch_ts), `TS`, `TSGetApplicationContext()`
2025 @*/
2026 PetscErrorCode TSSetApplicationContext(TS ts, void *ctx)
2027 {
2028   PetscFunctionBegin;
2029   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2030   ts->ctx = ctx;
2031   PetscFunctionReturn(PETSC_SUCCESS);
2032 }
2033 
2034 /*@
2035   TSGetApplicationContext - Gets the user-defined context for the
2036   timestepper that was set with `TSSetApplicationContext()`
2037 
2038   Not Collective
2039 
2040   Input Parameter:
2041 . ts - the `TS` context obtained from `TSCreate()`
2042 
2043   Output Parameter:
2044 . ctx - user context
2045 
2046   Level: intermediate
2047 
2048   Fortran Notes:
2049   You must write a Fortran interface definition for this
2050   function that tells Fortran the Fortran derived data type that you are passing in as the `ctx` argument.
2051 
2052 .seealso: [](ch_ts), `TS`, `TSSetApplicationContext()`
2053 @*/
2054 PetscErrorCode TSGetApplicationContext(TS ts, void *ctx)
2055 {
2056   PetscFunctionBegin;
2057   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2058   *(void **)ctx = ts->ctx;
2059   PetscFunctionReturn(PETSC_SUCCESS);
2060 }
2061 
2062 /*@
2063   TSGetStepNumber - Gets the number of time steps completed.
2064 
2065   Not Collective
2066 
2067   Input Parameter:
2068 . ts - the `TS` context obtained from `TSCreate()`
2069 
2070   Output Parameter:
2071 . steps - number of steps completed so far
2072 
2073   Level: intermediate
2074 
2075 .seealso: [](ch_ts), `TS`, `TSGetTime()`, `TSGetTimeStep()`, `TSSetPreStep()`, `TSSetPreStage()`, `TSSetPostStage()`, `TSSetPostStep()`
2076 @*/
2077 PetscErrorCode TSGetStepNumber(TS ts, PetscInt *steps)
2078 {
2079   PetscFunctionBegin;
2080   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2081   PetscAssertPointer(steps, 2);
2082   *steps = ts->steps;
2083   PetscFunctionReturn(PETSC_SUCCESS);
2084 }
2085 
2086 /*@
2087   TSSetStepNumber - Sets the number of steps completed.
2088 
2089   Logically Collective
2090 
2091   Input Parameters:
2092 + ts    - the `TS` context
2093 - steps - number of steps completed so far
2094 
2095   Level: developer
2096 
2097   Note:
2098   For most uses of the `TS` solvers the user need not explicitly call
2099   `TSSetStepNumber()`, as the step counter is appropriately updated in
2100   `TSSolve()`/`TSStep()`/`TSRollBack()`. Power users may call this routine to
2101   reinitialize timestepping by setting the step counter to zero (and time
2102   to the initial time) to solve a similar problem with different initial
2103   conditions or parameters. Other possible use case is to continue
2104   timestepping from a previously interrupted run in such a way that `TS`
2105   monitors will be called with a initial nonzero step counter.
2106 
2107 .seealso: [](ch_ts), `TS`, `TSGetStepNumber()`, `TSSetTime()`, `TSSetTimeStep()`, `TSSetSolution()`
2108 @*/
2109 PetscErrorCode TSSetStepNumber(TS ts, PetscInt steps)
2110 {
2111   PetscFunctionBegin;
2112   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2113   PetscValidLogicalCollectiveInt(ts, steps, 2);
2114   PetscCheck(steps >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Step number must be non-negative");
2115   ts->steps = steps;
2116   PetscFunctionReturn(PETSC_SUCCESS);
2117 }
2118 
2119 /*@
2120   TSSetTimeStep - Allows one to reset the timestep at any time,
2121   useful for simple pseudo-timestepping codes.
2122 
2123   Logically Collective
2124 
2125   Input Parameters:
2126 + ts        - the `TS` context obtained from `TSCreate()`
2127 - time_step - the size of the timestep
2128 
2129   Level: intermediate
2130 
2131 .seealso: [](ch_ts), `TS`, `TSPSEUDO`, `TSGetTimeStep()`, `TSSetTime()`
2132 @*/
2133 PetscErrorCode TSSetTimeStep(TS ts, PetscReal time_step)
2134 {
2135   PetscFunctionBegin;
2136   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2137   PetscValidLogicalCollectiveReal(ts, time_step, 2);
2138   ts->time_step = time_step;
2139   PetscFunctionReturn(PETSC_SUCCESS);
2140 }
2141 
2142 /*@
2143   TSSetExactFinalTime - Determines whether to adapt the final time step to
2144   match the exact final time, to interpolate the solution to the exact final time,
2145   or to just return at the final time `TS` computed (which may be slightly larger
2146   than the requested final time).
2147 
2148   Logically Collective
2149 
2150   Input Parameters:
2151 + ts     - the time-step context
2152 - eftopt - exact final time option
2153 .vb
2154   TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded, just use it
2155   TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time if the final time is exceeded
2156   TS_EXACTFINALTIME_MATCHSTEP   - Adapt final time step to ensure the computed final time exactly equals the requested final time
2157 .ve
2158 
2159   Options Database Key:
2160 . -ts_exact_final_time <stepover,interpolate,matchstep> - select the final step approach at runtime
2161 
2162   Level: beginner
2163 
2164   Note:
2165   If you use the option `TS_EXACTFINALTIME_STEPOVER` the solution may be at a very different time
2166   then the final time you selected.
2167 
2168 .seealso: [](ch_ts), `TS`, `TSExactFinalTimeOption`, `TSGetExactFinalTime()`
2169 @*/
2170 PetscErrorCode TSSetExactFinalTime(TS ts, TSExactFinalTimeOption eftopt)
2171 {
2172   PetscFunctionBegin;
2173   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2174   PetscValidLogicalCollectiveEnum(ts, eftopt, 2);
2175   ts->exact_final_time = eftopt;
2176   PetscFunctionReturn(PETSC_SUCCESS);
2177 }
2178 
2179 /*@
2180   TSGetExactFinalTime - Gets the exact final time option set with `TSSetExactFinalTime()`
2181 
2182   Not Collective
2183 
2184   Input Parameter:
2185 . ts - the `TS` context
2186 
2187   Output Parameter:
2188 . eftopt - exact final time option
2189 
2190   Level: beginner
2191 
2192 .seealso: [](ch_ts), `TS`, `TSExactFinalTimeOption`, `TSSetExactFinalTime()`
2193 @*/
2194 PetscErrorCode TSGetExactFinalTime(TS ts, TSExactFinalTimeOption *eftopt)
2195 {
2196   PetscFunctionBegin;
2197   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2198   PetscAssertPointer(eftopt, 2);
2199   *eftopt = ts->exact_final_time;
2200   PetscFunctionReturn(PETSC_SUCCESS);
2201 }
2202 
2203 /*@
2204   TSGetTimeStep - Gets the current timestep size.
2205 
2206   Not Collective
2207 
2208   Input Parameter:
2209 . ts - the `TS` context obtained from `TSCreate()`
2210 
2211   Output Parameter:
2212 . dt - the current timestep size
2213 
2214   Level: intermediate
2215 
2216 .seealso: [](ch_ts), `TS`, `TSSetTimeStep()`, `TSGetTime()`
2217 @*/
2218 PetscErrorCode TSGetTimeStep(TS ts, PetscReal *dt)
2219 {
2220   PetscFunctionBegin;
2221   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2222   PetscAssertPointer(dt, 2);
2223   *dt = ts->time_step;
2224   PetscFunctionReturn(PETSC_SUCCESS);
2225 }
2226 
2227 /*@
2228   TSGetSolution - Returns the solution at the present timestep. It
2229   is valid to call this routine inside the function that you are evaluating
2230   in order to move to the new timestep. This vector not changed until
2231   the solution at the next timestep has been calculated.
2232 
2233   Not Collective, but v returned is parallel if ts is parallel
2234 
2235   Input Parameter:
2236 . ts - the `TS` context obtained from `TSCreate()`
2237 
2238   Output Parameter:
2239 . v - the vector containing the solution
2240 
2241   Level: intermediate
2242 
2243   Note:
2244   If you used `TSSetExactFinalTime`(ts,`TS_EXACTFINALTIME_MATCHSTEP`); this does not return the solution at the requested
2245   final time. It returns the solution at the next timestep.
2246 
2247 .seealso: [](ch_ts), `TS`, `TSGetTimeStep()`, `TSGetTime()`, `TSGetSolveTime()`, `TSGetSolutionComponents()`, `TSSetSolutionFunction()`
2248 @*/
2249 PetscErrorCode TSGetSolution(TS ts, Vec *v)
2250 {
2251   PetscFunctionBegin;
2252   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2253   PetscAssertPointer(v, 2);
2254   *v = ts->vec_sol;
2255   PetscFunctionReturn(PETSC_SUCCESS);
2256 }
2257 
2258 /*@
2259   TSGetSolutionComponents - Returns any solution components at the present
2260   timestep, if available for the time integration method being used.
2261   Solution components are quantities that share the same size and
2262   structure as the solution vector.
2263 
2264   Not Collective, but v returned is parallel if ts is parallel
2265 
2266   Input Parameters:
2267 + ts - the `TS` context obtained from `TSCreate()` (input parameter).
2268 . n  - If v is `NULL`, then the number of solution components is
2269        returned through n, else the n-th solution component is
2270        returned in v.
2271 - v  - the vector containing the n-th solution component
2272        (may be `NULL` to use this function to find out
2273         the number of solutions components).
2274 
2275   Level: advanced
2276 
2277 .seealso: [](ch_ts), `TS`, `TSGetSolution()`
2278 @*/
2279 PetscErrorCode TSGetSolutionComponents(TS ts, PetscInt *n, Vec *v)
2280 {
2281   PetscFunctionBegin;
2282   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2283   if (!ts->ops->getsolutioncomponents) *n = 0;
2284   else PetscUseTypeMethod(ts, getsolutioncomponents, n, v);
2285   PetscFunctionReturn(PETSC_SUCCESS);
2286 }
2287 
2288 /*@
2289   TSGetAuxSolution - Returns an auxiliary solution at the present
2290   timestep, if available for the time integration method being used.
2291 
2292   Not Collective, but v returned is parallel if ts is parallel
2293 
2294   Input Parameters:
2295 + ts - the `TS` context obtained from `TSCreate()` (input parameter).
2296 - v  - the vector containing the auxiliary solution
2297 
2298   Level: intermediate
2299 
2300 .seealso: [](ch_ts), `TS`, `TSGetSolution()`
2301 @*/
2302 PetscErrorCode TSGetAuxSolution(TS ts, Vec *v)
2303 {
2304   PetscFunctionBegin;
2305   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2306   if (ts->ops->getauxsolution) PetscUseTypeMethod(ts, getauxsolution, v);
2307   else PetscCall(VecZeroEntries(*v));
2308   PetscFunctionReturn(PETSC_SUCCESS);
2309 }
2310 
2311 /*@
2312   TSGetTimeError - Returns the estimated error vector, if the chosen
2313   `TSType` has an error estimation functionality and `TSSetTimeError()` was called
2314 
2315   Not Collective, but v returned is parallel if ts is parallel
2316 
2317   Input Parameters:
2318 + ts - the `TS` context obtained from `TSCreate()` (input parameter).
2319 . n  - current estimate (n=0) or previous one (n=-1)
2320 - v  - the vector containing the error (same size as the solution).
2321 
2322   Level: intermediate
2323 
2324   Note:
2325   MUST call after `TSSetUp()`
2326 
2327 .seealso: [](ch_ts), `TSGetSolution()`, `TSSetTimeError()`
2328 @*/
2329 PetscErrorCode TSGetTimeError(TS ts, PetscInt n, Vec *v)
2330 {
2331   PetscFunctionBegin;
2332   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2333   if (ts->ops->gettimeerror) PetscUseTypeMethod(ts, gettimeerror, n, v);
2334   else PetscCall(VecZeroEntries(*v));
2335   PetscFunctionReturn(PETSC_SUCCESS);
2336 }
2337 
2338 /*@
2339   TSSetTimeError - Sets the estimated error vector, if the chosen
2340   `TSType` has an error estimation functionality. This can be used
2341   to restart such a time integrator with a given error vector.
2342 
2343   Not Collective, but v returned is parallel if ts is parallel
2344 
2345   Input Parameters:
2346 + ts - the `TS` context obtained from `TSCreate()` (input parameter).
2347 - v  - the vector containing the error (same size as the solution).
2348 
2349   Level: intermediate
2350 
2351 .seealso: [](ch_ts), `TS`, `TSSetSolution()`, `TSGetTimeError()`
2352 @*/
2353 PetscErrorCode TSSetTimeError(TS ts, Vec v)
2354 {
2355   PetscFunctionBegin;
2356   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2357   PetscCheck(ts->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetUp() first");
2358   PetscTryTypeMethod(ts, settimeerror, v);
2359   PetscFunctionReturn(PETSC_SUCCESS);
2360 }
2361 
2362 /* ----- Routines to initialize and destroy a timestepper ---- */
2363 /*@
2364   TSSetProblemType - Sets the type of problem to be solved.
2365 
2366   Not collective
2367 
2368   Input Parameters:
2369 + ts   - The `TS`
2370 - type - One of `TS_LINEAR`, `TS_NONLINEAR` where these types refer to problems of the forms
2371 .vb
2372          U_t - A U = 0      (linear)
2373          U_t - A(t) U = 0   (linear)
2374          F(t,U,U_t) = 0     (nonlinear)
2375 .ve
2376 
2377   Level: beginner
2378 
2379 .seealso: [](ch_ts), `TSSetUp()`, `TSProblemType`, `TS`
2380 @*/
2381 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type)
2382 {
2383   PetscFunctionBegin;
2384   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2385   ts->problem_type = type;
2386   if (type == TS_LINEAR) {
2387     SNES snes;
2388     PetscCall(TSGetSNES(ts, &snes));
2389     PetscCall(SNESSetType(snes, SNESKSPONLY));
2390   }
2391   PetscFunctionReturn(PETSC_SUCCESS);
2392 }
2393 
2394 /*@
2395   TSGetProblemType - Gets the type of problem to be solved.
2396 
2397   Not collective
2398 
2399   Input Parameter:
2400 . ts - The `TS`
2401 
2402   Output Parameter:
2403 . type - One of `TS_LINEAR`, `TS_NONLINEAR` where these types refer to problems of the forms
2404 .vb
2405          M U_t = A U
2406          M(t) U_t = A(t) U
2407          F(t,U,U_t)
2408 .ve
2409 
2410   Level: beginner
2411 
2412 .seealso: [](ch_ts), `TSSetUp()`, `TSProblemType`, `TS`
2413 @*/
2414 PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type)
2415 {
2416   PetscFunctionBegin;
2417   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2418   PetscAssertPointer(type, 2);
2419   *type = ts->problem_type;
2420   PetscFunctionReturn(PETSC_SUCCESS);
2421 }
2422 
2423 /*
2424     Attempt to check/preset a default value for the exact final time option. This is needed at the beginning of TSSolve() and in TSSetUp()
2425 */
2426 static PetscErrorCode TSSetExactFinalTimeDefault(TS ts)
2427 {
2428   PetscBool isnone;
2429 
2430   PetscFunctionBegin;
2431   PetscCall(TSGetAdapt(ts, &ts->adapt));
2432   PetscCall(TSAdaptSetDefaultType(ts->adapt, ts->default_adapt_type));
2433 
2434   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &isnone));
2435   if (!isnone && ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
2436   else if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) ts->exact_final_time = TS_EXACTFINALTIME_INTERPOLATE;
2437   PetscFunctionReturn(PETSC_SUCCESS);
2438 }
2439 
2440 /*@
2441   TSSetUp - Sets up the internal data structures for the later use of a timestepper.
2442 
2443   Collective
2444 
2445   Input Parameter:
2446 . ts - the `TS` context obtained from `TSCreate()`
2447 
2448   Level: advanced
2449 
2450   Note:
2451   For basic use of the `TS` solvers the user need not explicitly call
2452   `TSSetUp()`, since these actions will automatically occur during
2453   the call to `TSStep()` or `TSSolve()`.  However, if one wishes to control this
2454   phase separately, `TSSetUp()` should be called after `TSCreate()`
2455   and optional routines of the form TSSetXXX(), but before `TSStep()` and `TSSolve()`.
2456 
2457 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSStep()`, `TSDestroy()`, `TSSolve()`
2458 @*/
2459 PetscErrorCode TSSetUp(TS ts)
2460 {
2461   DM dm;
2462   PetscErrorCode (*func)(SNES, Vec, Vec, void *);
2463   PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *);
2464   TSIFunctionFn   *ifun;
2465   TSIJacobianFn   *ijac;
2466   TSI2JacobianFn  *i2jac;
2467   TSRHSJacobianFn *rhsjac;
2468 
2469   PetscFunctionBegin;
2470   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2471   if (ts->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
2472 
2473   if (!((PetscObject)ts)->type_name) {
2474     PetscCall(TSGetIFunction(ts, NULL, &ifun, NULL));
2475     PetscCall(TSSetType(ts, ifun ? TSBEULER : TSEULER));
2476   }
2477 
2478   if (!ts->vec_sol) {
2479     PetscCheck(ts->dm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetSolution() first");
2480     PetscCall(DMCreateGlobalVector(ts->dm, &ts->vec_sol));
2481   }
2482 
2483   if (ts->eval_times) {
2484     if (!ts->eval_times->sol_vecs) PetscCall(VecDuplicateVecs(ts->vec_sol, ts->eval_times->num_time_points, &ts->eval_times->sol_vecs));
2485     if (!ts->eval_times->sol_times) PetscCall(PetscMalloc1(ts->eval_times->num_time_points, &ts->eval_times->sol_times));
2486   }
2487   if (!ts->Jacp && ts->Jacprhs) { /* IJacobianP shares the same matrix with RHSJacobianP if only RHSJacobianP is provided */
2488     PetscCall(PetscObjectReference((PetscObject)ts->Jacprhs));
2489     ts->Jacp = ts->Jacprhs;
2490   }
2491 
2492   if (ts->quadraturets) {
2493     PetscCall(TSSetUp(ts->quadraturets));
2494     PetscCall(VecDestroy(&ts->vec_costintegrand));
2495     PetscCall(VecDuplicate(ts->quadraturets->vec_sol, &ts->vec_costintegrand));
2496   }
2497 
2498   PetscCall(TSGetRHSJacobian(ts, NULL, NULL, &rhsjac, NULL));
2499   if (rhsjac == TSComputeRHSJacobianConstant) {
2500     Mat  Amat, Pmat;
2501     SNES snes;
2502     PetscCall(TSGetSNES(ts, &snes));
2503     PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL));
2504     /* Matching matrices implies that an IJacobian is NOT set, because if it had been set, the IJacobian's matrix would
2505      * have displaced the RHS matrix */
2506     if (Amat && Amat == ts->Arhs) {
2507       /* we need to copy the values of the matrix because for the constant Jacobian case the user will never set the numerical values in this new location */
2508       PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
2509       PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL));
2510       PetscCall(MatDestroy(&Amat));
2511     }
2512     if (Pmat && Pmat == ts->Brhs) {
2513       PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat));
2514       PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL));
2515       PetscCall(MatDestroy(&Pmat));
2516     }
2517   }
2518 
2519   PetscCall(TSGetAdapt(ts, &ts->adapt));
2520   PetscCall(TSAdaptSetDefaultType(ts->adapt, ts->default_adapt_type));
2521 
2522   PetscTryTypeMethod(ts, setup);
2523 
2524   PetscCall(TSSetExactFinalTimeDefault(ts));
2525 
2526   /* In the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
2527      to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
2528    */
2529   PetscCall(TSGetDM(ts, &dm));
2530   PetscCall(DMSNESGetFunction(dm, &func, NULL));
2531   if (!func) PetscCall(DMSNESSetFunction(dm, SNESTSFormFunction, ts));
2532 
2533   /* If the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
2534      Otherwise, the SNES will use coloring internally to form the Jacobian.
2535    */
2536   PetscCall(DMSNESGetJacobian(dm, &jac, NULL));
2537   PetscCall(DMTSGetIJacobian(dm, &ijac, NULL));
2538   PetscCall(DMTSGetI2Jacobian(dm, &i2jac, NULL));
2539   PetscCall(DMTSGetRHSJacobian(dm, &rhsjac, NULL));
2540   if (!jac && (ijac || i2jac || rhsjac)) PetscCall(DMSNESSetJacobian(dm, SNESTSFormJacobian, ts));
2541 
2542   /* if time integration scheme has a starting method, call it */
2543   PetscTryTypeMethod(ts, startingmethod);
2544 
2545   ts->setupcalled = PETSC_TRUE;
2546   PetscFunctionReturn(PETSC_SUCCESS);
2547 }
2548 
2549 /*@
2550   TSReset - Resets a `TS` context to the state it was in before `TSSetUp()` was called and removes any allocated `Vec` and `Mat` from its data structures
2551 
2552   Collective
2553 
2554   Input Parameter:
2555 . ts - the `TS` context obtained from `TSCreate()`
2556 
2557   Level: developer
2558 
2559   Notes:
2560   Any options set on the `TS` object, including those set with `TSSetFromOptions()` remain.
2561 
2562   See also `TSSetResize()` to change the size of the system being integrated (for example by adaptive mesh refinement) during the time integration.
2563 
2564 .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetup()`, `TSDestroy()`, `TSSetResize()`
2565 @*/
2566 PetscErrorCode TSReset(TS ts)
2567 {
2568   TS_RHSSplitLink ilink = ts->tsrhssplit, next;
2569 
2570   PetscFunctionBegin;
2571   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2572 
2573   PetscTryTypeMethod(ts, reset);
2574   if (ts->snes) PetscCall(SNESReset(ts->snes));
2575   if (ts->adapt) PetscCall(TSAdaptReset(ts->adapt));
2576 
2577   PetscCall(MatDestroy(&ts->Arhs));
2578   PetscCall(MatDestroy(&ts->Brhs));
2579   PetscCall(VecDestroy(&ts->Frhs));
2580   PetscCall(VecDestroy(&ts->vec_sol));
2581   PetscCall(VecDestroy(&ts->vec_sol0));
2582   PetscCall(VecDestroy(&ts->vec_dot));
2583   PetscCall(VecDestroy(&ts->vatol));
2584   PetscCall(VecDestroy(&ts->vrtol));
2585   PetscCall(VecDestroyVecs(ts->nwork, &ts->work));
2586 
2587   PetscCall(MatDestroy(&ts->Jacprhs));
2588   PetscCall(MatDestroy(&ts->Jacp));
2589   if (ts->forward_solve) PetscCall(TSForwardReset(ts));
2590   if (ts->quadraturets) {
2591     PetscCall(TSReset(ts->quadraturets));
2592     PetscCall(VecDestroy(&ts->vec_costintegrand));
2593   }
2594   while (ilink) {
2595     next = ilink->next;
2596     PetscCall(TSDestroy(&ilink->ts));
2597     PetscCall(PetscFree(ilink->splitname));
2598     PetscCall(ISDestroy(&ilink->is));
2599     PetscCall(PetscFree(ilink));
2600     ilink = next;
2601   }
2602   ts->tsrhssplit     = NULL;
2603   ts->num_rhs_splits = 0;
2604   if (ts->eval_times) {
2605     PetscCall(PetscFree(ts->eval_times->time_points));
2606     PetscCall(PetscFree(ts->eval_times->sol_times));
2607     PetscCall(VecDestroyVecs(ts->eval_times->num_time_points, &ts->eval_times->sol_vecs));
2608     PetscCall(PetscFree(ts->eval_times));
2609   }
2610   ts->rhsjacobian.time  = PETSC_MIN_REAL;
2611   ts->rhsjacobian.scale = 1.0;
2612   ts->ijacobian.shift   = 1.0;
2613   ts->setupcalled       = PETSC_FALSE;
2614   PetscFunctionReturn(PETSC_SUCCESS);
2615 }
2616 
2617 static PetscErrorCode TSResizeReset(TS);
2618 
2619 /*@
2620   TSDestroy - Destroys the timestepper context that was created
2621   with `TSCreate()`.
2622 
2623   Collective
2624 
2625   Input Parameter:
2626 . ts - the `TS` context obtained from `TSCreate()`
2627 
2628   Level: beginner
2629 
2630 .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetUp()`, `TSSolve()`
2631 @*/
2632 PetscErrorCode TSDestroy(TS *ts)
2633 {
2634   PetscFunctionBegin;
2635   if (!*ts) PetscFunctionReturn(PETSC_SUCCESS);
2636   PetscValidHeaderSpecific(*ts, TS_CLASSID, 1);
2637   if (--((PetscObject)*ts)->refct > 0) {
2638     *ts = NULL;
2639     PetscFunctionReturn(PETSC_SUCCESS);
2640   }
2641 
2642   PetscCall(TSReset(*ts));
2643   PetscCall(TSAdjointReset(*ts));
2644   if ((*ts)->forward_solve) PetscCall(TSForwardReset(*ts));
2645   PetscCall(TSResizeReset(*ts));
2646 
2647   /* if memory was published with SAWs then destroy it */
2648   PetscCall(PetscObjectSAWsViewOff((PetscObject)*ts));
2649   PetscTryTypeMethod(*ts, destroy);
2650 
2651   PetscCall(TSTrajectoryDestroy(&(*ts)->trajectory));
2652 
2653   PetscCall(TSAdaptDestroy(&(*ts)->adapt));
2654   PetscCall(TSEventDestroy(&(*ts)->event));
2655 
2656   PetscCall(SNESDestroy(&(*ts)->snes));
2657   PetscCall(SNESDestroy(&(*ts)->snesrhssplit));
2658   PetscCall(DMDestroy(&(*ts)->dm));
2659   PetscCall(TSMonitorCancel(*ts));
2660   PetscCall(TSAdjointMonitorCancel(*ts));
2661 
2662   PetscCall(TSDestroy(&(*ts)->quadraturets));
2663   PetscCall(PetscHeaderDestroy(ts));
2664   PetscFunctionReturn(PETSC_SUCCESS);
2665 }
2666 
2667 /*@
2668   TSGetSNES - Returns the `SNES` (nonlinear solver) associated with
2669   a `TS` (timestepper) context. Valid only for nonlinear problems.
2670 
2671   Not Collective, but snes is parallel if ts is parallel
2672 
2673   Input Parameter:
2674 . ts - the `TS` context obtained from `TSCreate()`
2675 
2676   Output Parameter:
2677 . snes - the nonlinear solver context
2678 
2679   Level: beginner
2680 
2681   Notes:
2682   The user can then directly manipulate the `SNES` context to set various
2683   options, etc.  Likewise, the user can then extract and manipulate the
2684   `KSP`, and `PC` contexts as well.
2685 
2686   `TSGetSNES()` does not work for integrators that do not use `SNES`; in
2687   this case `TSGetSNES()` returns `NULL` in `snes`.
2688 
2689 .seealso: [](ch_ts), `TS`, `SNES`, `TSCreate()`, `TSSetUp()`, `TSSolve()`
2690 @*/
2691 PetscErrorCode TSGetSNES(TS ts, SNES *snes)
2692 {
2693   PetscFunctionBegin;
2694   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2695   PetscAssertPointer(snes, 2);
2696   if (!ts->snes) {
2697     PetscCall(SNESCreate(PetscObjectComm((PetscObject)ts), &ts->snes));
2698     PetscCall(PetscObjectSetOptions((PetscObject)ts->snes, ((PetscObject)ts)->options));
2699     PetscCall(SNESSetFunction(ts->snes, NULL, SNESTSFormFunction, ts));
2700     PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->snes, (PetscObject)ts, 1));
2701     if (ts->dm) PetscCall(SNESSetDM(ts->snes, ts->dm));
2702     if (ts->problem_type == TS_LINEAR) PetscCall(SNESSetType(ts->snes, SNESKSPONLY));
2703   }
2704   *snes = ts->snes;
2705   PetscFunctionReturn(PETSC_SUCCESS);
2706 }
2707 
2708 /*@
2709   TSSetSNES - Set the `SNES` (nonlinear solver) to be used by the `TS` timestepping context
2710 
2711   Collective
2712 
2713   Input Parameters:
2714 + ts   - the `TS` context obtained from `TSCreate()`
2715 - snes - the nonlinear solver context
2716 
2717   Level: developer
2718 
2719   Note:
2720   Most users should have the `TS` created by calling `TSGetSNES()`
2721 
2722 .seealso: [](ch_ts), `TS`, `SNES`, `TSCreate()`, `TSSetUp()`, `TSSolve()`, `TSGetSNES()`
2723 @*/
2724 PetscErrorCode TSSetSNES(TS ts, SNES snes)
2725 {
2726   PetscErrorCode (*func)(SNES, Vec, Mat, Mat, void *);
2727 
2728   PetscFunctionBegin;
2729   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2730   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
2731   PetscCall(PetscObjectReference((PetscObject)snes));
2732   PetscCall(SNESDestroy(&ts->snes));
2733 
2734   ts->snes = snes;
2735 
2736   PetscCall(SNESSetFunction(ts->snes, NULL, SNESTSFormFunction, ts));
2737   PetscCall(SNESGetJacobian(ts->snes, NULL, NULL, &func, NULL));
2738   if (func == SNESTSFormJacobian) PetscCall(SNESSetJacobian(ts->snes, NULL, NULL, SNESTSFormJacobian, ts));
2739   PetscFunctionReturn(PETSC_SUCCESS);
2740 }
2741 
2742 /*@
2743   TSGetKSP - Returns the `KSP` (linear solver) associated with
2744   a `TS` (timestepper) context.
2745 
2746   Not Collective, but `ksp` is parallel if `ts` is parallel
2747 
2748   Input Parameter:
2749 . ts - the `TS` context obtained from `TSCreate()`
2750 
2751   Output Parameter:
2752 . ksp - the nonlinear solver context
2753 
2754   Level: beginner
2755 
2756   Notes:
2757   The user can then directly manipulate the `KSP` context to set various
2758   options, etc.  Likewise, the user can then extract and manipulate the
2759   `PC` context as well.
2760 
2761   `TSGetKSP()` does not work for integrators that do not use `KSP`;
2762   in this case `TSGetKSP()` returns `NULL` in `ksp`.
2763 
2764 .seealso: [](ch_ts), `TS`, `SNES`, `KSP`, `TSCreate()`, `TSSetUp()`, `TSSolve()`, `TSGetSNES()`
2765 @*/
2766 PetscErrorCode TSGetKSP(TS ts, KSP *ksp)
2767 {
2768   SNES snes;
2769 
2770   PetscFunctionBegin;
2771   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2772   PetscAssertPointer(ksp, 2);
2773   PetscCheck(((PetscObject)ts)->type_name, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "KSP is not created yet. Call TSSetType() first");
2774   PetscCheck(ts->problem_type == TS_LINEAR, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Linear only; use TSGetSNES()");
2775   PetscCall(TSGetSNES(ts, &snes));
2776   PetscCall(SNESGetKSP(snes, ksp));
2777   PetscFunctionReturn(PETSC_SUCCESS);
2778 }
2779 
2780 /* ----------- Routines to set solver parameters ---------- */
2781 
2782 /*@
2783   TSSetMaxSteps - Sets the maximum number of steps to use.
2784 
2785   Logically Collective
2786 
2787   Input Parameters:
2788 + ts       - the `TS` context obtained from `TSCreate()`
2789 - maxsteps - maximum number of steps to use
2790 
2791   Options Database Key:
2792 . -ts_max_steps <maxsteps> - Sets maxsteps
2793 
2794   Level: intermediate
2795 
2796   Note:
2797   Use `PETSC_DETERMINE` to reset the maximum number of steps to the default from when the object's type was set
2798 
2799   The default maximum number of steps is 5,000
2800 
2801   Fortran Note:
2802   Use `PETSC_DETERMINE_INTEGER`
2803 
2804 .seealso: [](ch_ts), `TS`, `TSGetMaxSteps()`, `TSSetMaxTime()`, `TSSetExactFinalTime()`
2805 @*/
2806 PetscErrorCode TSSetMaxSteps(TS ts, PetscInt maxsteps)
2807 {
2808   PetscFunctionBegin;
2809   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2810   PetscValidLogicalCollectiveInt(ts, maxsteps, 2);
2811   if (maxsteps == PETSC_DETERMINE) {
2812     ts->max_steps = ts->default_max_steps;
2813   } else {
2814     PetscCheck(maxsteps >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of steps must be non-negative");
2815     ts->max_steps = maxsteps;
2816   }
2817   PetscFunctionReturn(PETSC_SUCCESS);
2818 }
2819 
2820 /*@
2821   TSGetMaxSteps - Gets the maximum number of steps to use.
2822 
2823   Not Collective
2824 
2825   Input Parameter:
2826 . ts - the `TS` context obtained from `TSCreate()`
2827 
2828   Output Parameter:
2829 . maxsteps - maximum number of steps to use
2830 
2831   Level: advanced
2832 
2833 .seealso: [](ch_ts), `TS`, `TSSetMaxSteps()`, `TSGetMaxTime()`, `TSSetMaxTime()`
2834 @*/
2835 PetscErrorCode TSGetMaxSteps(TS ts, PetscInt *maxsteps)
2836 {
2837   PetscFunctionBegin;
2838   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2839   PetscAssertPointer(maxsteps, 2);
2840   *maxsteps = ts->max_steps;
2841   PetscFunctionReturn(PETSC_SUCCESS);
2842 }
2843 
2844 /*@
2845   TSSetMaxTime - Sets the maximum (or final) time for timestepping.
2846 
2847   Logically Collective
2848 
2849   Input Parameters:
2850 + ts      - the `TS` context obtained from `TSCreate()`
2851 - maxtime - final time to step to
2852 
2853   Options Database Key:
2854 . -ts_max_time <maxtime> - Sets maxtime
2855 
2856   Level: intermediate
2857 
2858   Notes:
2859   Use `PETSC_DETERMINE` to reset the maximum time to the default from when the object's type was set
2860 
2861   The default maximum time is 5.0
2862 
2863   Fortran Note:
2864   Use `PETSC_DETERMINE_REAL`
2865 
2866 .seealso: [](ch_ts), `TS`, `TSGetMaxTime()`, `TSSetMaxSteps()`, `TSSetExactFinalTime()`
2867 @*/
2868 PetscErrorCode TSSetMaxTime(TS ts, PetscReal maxtime)
2869 {
2870   PetscFunctionBegin;
2871   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2872   PetscValidLogicalCollectiveReal(ts, maxtime, 2);
2873   if (maxtime == PETSC_DETERMINE) {
2874     ts->max_time = ts->default_max_time;
2875   } else {
2876     ts->max_time = maxtime;
2877   }
2878   PetscFunctionReturn(PETSC_SUCCESS);
2879 }
2880 
2881 /*@
2882   TSGetMaxTime - Gets the maximum (or final) time for timestepping.
2883 
2884   Not Collective
2885 
2886   Input Parameter:
2887 . ts - the `TS` context obtained from `TSCreate()`
2888 
2889   Output Parameter:
2890 . maxtime - final time to step to
2891 
2892   Level: advanced
2893 
2894 .seealso: [](ch_ts), `TS`, `TSSetMaxTime()`, `TSGetMaxSteps()`, `TSSetMaxSteps()`
2895 @*/
2896 PetscErrorCode TSGetMaxTime(TS ts, PetscReal *maxtime)
2897 {
2898   PetscFunctionBegin;
2899   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2900   PetscAssertPointer(maxtime, 2);
2901   *maxtime = ts->max_time;
2902   PetscFunctionReturn(PETSC_SUCCESS);
2903 }
2904 
2905 // PetscClangLinter pragma disable: -fdoc-*
2906 /*@
2907   TSSetInitialTimeStep - Deprecated, use `TSSetTime()` and `TSSetTimeStep()`.
2908 
2909   Level: deprecated
2910 
2911 @*/
2912 PetscErrorCode TSSetInitialTimeStep(TS ts, PetscReal initial_time, PetscReal time_step)
2913 {
2914   PetscFunctionBegin;
2915   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2916   PetscCall(TSSetTime(ts, initial_time));
2917   PetscCall(TSSetTimeStep(ts, time_step));
2918   PetscFunctionReturn(PETSC_SUCCESS);
2919 }
2920 
2921 // PetscClangLinter pragma disable: -fdoc-*
2922 /*@
2923   TSGetDuration - Deprecated, use `TSGetMaxSteps()` and `TSGetMaxTime()`.
2924 
2925   Level: deprecated
2926 
2927 @*/
2928 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
2929 {
2930   PetscFunctionBegin;
2931   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2932   if (maxsteps) {
2933     PetscAssertPointer(maxsteps, 2);
2934     *maxsteps = ts->max_steps;
2935   }
2936   if (maxtime) {
2937     PetscAssertPointer(maxtime, 3);
2938     *maxtime = ts->max_time;
2939   }
2940   PetscFunctionReturn(PETSC_SUCCESS);
2941 }
2942 
2943 // PetscClangLinter pragma disable: -fdoc-*
2944 /*@
2945   TSSetDuration - Deprecated, use `TSSetMaxSteps()` and `TSSetMaxTime()`.
2946 
2947   Level: deprecated
2948 
2949 @*/
2950 PetscErrorCode TSSetDuration(TS ts, PetscInt maxsteps, PetscReal maxtime)
2951 {
2952   PetscFunctionBegin;
2953   if (maxsteps != PETSC_CURRENT) PetscCall(TSSetMaxSteps(ts, maxsteps));
2954   if (maxtime != (PetscReal)PETSC_CURRENT) PetscCall(TSSetMaxTime(ts, maxtime));
2955   PetscFunctionReturn(PETSC_SUCCESS);
2956 }
2957 
2958 // PetscClangLinter pragma disable: -fdoc-*
2959 /*@
2960   TSGetTimeStepNumber - Deprecated, use `TSGetStepNumber()`.
2961 
2962   Level: deprecated
2963 
2964 @*/
2965 PetscErrorCode TSGetTimeStepNumber(TS ts, PetscInt *steps)
2966 {
2967   return TSGetStepNumber(ts, steps);
2968 }
2969 
2970 // PetscClangLinter pragma disable: -fdoc-*
2971 /*@
2972   TSGetTotalSteps - Deprecated, use `TSGetStepNumber()`.
2973 
2974   Level: deprecated
2975 
2976 @*/
2977 PetscErrorCode TSGetTotalSteps(TS ts, PetscInt *steps)
2978 {
2979   return TSGetStepNumber(ts, steps);
2980 }
2981 
2982 /*@
2983   TSSetSolution - Sets the initial solution vector
2984   for use by the `TS` routines.
2985 
2986   Logically Collective
2987 
2988   Input Parameters:
2989 + ts - the `TS` context obtained from `TSCreate()`
2990 - u  - the solution vector
2991 
2992   Level: beginner
2993 
2994 .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `TSGetSolution()`, `TSCreate()`
2995 @*/
2996 PetscErrorCode TSSetSolution(TS ts, Vec u)
2997 {
2998   DM dm;
2999 
3000   PetscFunctionBegin;
3001   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3002   PetscValidHeaderSpecific(u, VEC_CLASSID, 2);
3003   PetscCall(PetscObjectReference((PetscObject)u));
3004   PetscCall(VecDestroy(&ts->vec_sol));
3005   ts->vec_sol = u;
3006 
3007   PetscCall(TSGetDM(ts, &dm));
3008   PetscCall(DMShellSetGlobalVector(dm, u));
3009   PetscFunctionReturn(PETSC_SUCCESS);
3010 }
3011 
3012 /*@C
3013   TSSetPreStep - Sets the general-purpose function
3014   called once at the beginning of each time step.
3015 
3016   Logically Collective
3017 
3018   Input Parameters:
3019 + ts   - The `TS` context obtained from `TSCreate()`
3020 - func - The function
3021 
3022   Calling sequence of `func`:
3023 . ts - the `TS` context
3024 
3025   Level: intermediate
3026 
3027 .seealso: [](ch_ts), `TS`, `TSSetPreStage()`, `TSSetPostStage()`, `TSSetPostStep()`, `TSStep()`, `TSRestartStep()`
3028 @*/
3029 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS ts))
3030 {
3031   PetscFunctionBegin;
3032   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3033   ts->prestep = func;
3034   PetscFunctionReturn(PETSC_SUCCESS);
3035 }
3036 
3037 /*@
3038   TSPreStep - Runs the user-defined pre-step function provided with `TSSetPreStep()`
3039 
3040   Collective
3041 
3042   Input Parameter:
3043 . ts - The `TS` context obtained from `TSCreate()`
3044 
3045   Level: developer
3046 
3047   Note:
3048   `TSPreStep()` is typically used within time stepping implementations,
3049   so most users would not generally call this routine themselves.
3050 
3051 .seealso: [](ch_ts), `TS`, `TSSetPreStep()`, `TSPreStage()`, `TSPostStage()`, `TSPostStep()`
3052 @*/
3053 PetscErrorCode TSPreStep(TS ts)
3054 {
3055   PetscFunctionBegin;
3056   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3057   if (ts->prestep) {
3058     Vec              U;
3059     PetscObjectId    idprev;
3060     PetscBool        sameObject;
3061     PetscObjectState sprev, spost;
3062 
3063     PetscCall(TSGetSolution(ts, &U));
3064     PetscCall(PetscObjectGetId((PetscObject)U, &idprev));
3065     PetscCall(PetscObjectStateGet((PetscObject)U, &sprev));
3066     PetscCallBack("TS callback preset", (*ts->prestep)(ts));
3067     PetscCall(TSGetSolution(ts, &U));
3068     PetscCall(PetscObjectCompareId((PetscObject)U, idprev, &sameObject));
3069     PetscCall(PetscObjectStateGet((PetscObject)U, &spost));
3070     if (!sameObject || sprev != spost) PetscCall(TSRestartStep(ts));
3071   }
3072   PetscFunctionReturn(PETSC_SUCCESS);
3073 }
3074 
3075 /*@C
3076   TSSetPreStage - Sets the general-purpose function
3077   called once at the beginning of each stage.
3078 
3079   Logically Collective
3080 
3081   Input Parameters:
3082 + ts   - The `TS` context obtained from `TSCreate()`
3083 - func - The function
3084 
3085   Calling sequence of `func`:
3086 + ts        - the `TS` context
3087 - stagetime - the stage time
3088 
3089   Level: intermediate
3090 
3091   Note:
3092   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
3093   The time step number being computed can be queried using `TSGetStepNumber()` and the total size of the step being
3094   attempted can be obtained using `TSGetTimeStep()`. The time at the start of the step is available via `TSGetTime()`.
3095 
3096 .seealso: [](ch_ts), `TS`, `TSSetPostStage()`, `TSSetPreStep()`, `TSSetPostStep()`, `TSGetApplicationContext()`
3097 @*/
3098 PetscErrorCode TSSetPreStage(TS ts, PetscErrorCode (*func)(TS ts, PetscReal stagetime))
3099 {
3100   PetscFunctionBegin;
3101   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3102   ts->prestage = func;
3103   PetscFunctionReturn(PETSC_SUCCESS);
3104 }
3105 
3106 /*@C
3107   TSSetPostStage - Sets the general-purpose function
3108   called once at the end of each stage.
3109 
3110   Logically Collective
3111 
3112   Input Parameters:
3113 + ts   - The `TS` context obtained from `TSCreate()`
3114 - func - The function
3115 
3116   Calling sequence of `func`:
3117 + ts         - the `TS` context
3118 . stagetime  - the stage time
3119 . stageindex - the stage index
3120 - Y          - Array of vectors (of size = total number of stages) with the stage solutions
3121 
3122   Level: intermediate
3123 
3124   Note:
3125   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
3126   The time step number being computed can be queried using `TSGetStepNumber()` and the total size of the step being
3127   attempted can be obtained using `TSGetTimeStep()`. The time at the start of the step is available via `TSGetTime()`.
3128 
3129 .seealso: [](ch_ts), `TS`, `TSSetPreStage()`, `TSSetPreStep()`, `TSSetPostStep()`, `TSGetApplicationContext()`
3130 @*/
3131 PetscErrorCode TSSetPostStage(TS ts, PetscErrorCode (*func)(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y))
3132 {
3133   PetscFunctionBegin;
3134   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3135   ts->poststage = func;
3136   PetscFunctionReturn(PETSC_SUCCESS);
3137 }
3138 
3139 /*@C
3140   TSSetPostEvaluate - Sets the general-purpose function
3141   called at the end of each step evaluation.
3142 
3143   Logically Collective
3144 
3145   Input Parameters:
3146 + ts   - The `TS` context obtained from `TSCreate()`
3147 - func - The function
3148 
3149   Calling sequence of `func`:
3150 . ts - the `TS` context
3151 
3152   Level: intermediate
3153 
3154   Note:
3155   The function set by `TSSetPostEvaluate()` is called after the solution is evaluated, or after the step rollback.
3156   Inside the `func` callback, the solution vector can be obtained with `TSGetSolution()`, and modified, if need be.
3157   The time step can be obtained with `TSGetTimeStep()`, and the time at the start of the step - via `TSGetTime()`.
3158   The potential changes to the solution vector introduced by event handling (`postevent()`) are not relevant for `TSSetPostEvaluate()`,
3159   but are relevant for `TSSetPostStep()`, according to the function call scheme in `TSSolve()`, as shown below
3160 .vb
3161   ...
3162   Step()
3163   PostEvaluate()
3164   EventHandling()
3165   step_rollback ? PostEvaluate() : PostStep()
3166   ...
3167 .ve
3168   where EventHandling() may result in one of the following three outcomes
3169 .vb
3170   (1) | successful step | solution intact
3171   (2) | successful step | solution modified by `postevent()`
3172   (3) | step_rollback   | solution rolled back
3173 .ve
3174 
3175 .seealso: [](ch_ts), `TS`, `TSSetPreStage()`, `TSSetPreStep()`, `TSSetPostStep()`, `TSGetApplicationContext()`
3176 @*/
3177 PetscErrorCode TSSetPostEvaluate(TS ts, PetscErrorCode (*func)(TS ts))
3178 {
3179   PetscFunctionBegin;
3180   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3181   ts->postevaluate = func;
3182   PetscFunctionReturn(PETSC_SUCCESS);
3183 }
3184 
3185 /*@
3186   TSPreStage - Runs the user-defined pre-stage function set using `TSSetPreStage()`
3187 
3188   Collective
3189 
3190   Input Parameters:
3191 + ts        - The `TS` context obtained from `TSCreate()`
3192 - stagetime - The absolute time of the current stage
3193 
3194   Level: developer
3195 
3196   Note:
3197   `TSPreStage()` is typically used within time stepping implementations,
3198   most users would not generally call this routine themselves.
3199 
3200 .seealso: [](ch_ts), `TS`, `TSPostStage()`, `TSSetPreStep()`, `TSPreStep()`, `TSPostStep()`
3201 @*/
3202 PetscErrorCode TSPreStage(TS ts, PetscReal stagetime)
3203 {
3204   PetscFunctionBegin;
3205   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3206   if (ts->prestage) PetscCallBack("TS callback prestage", (*ts->prestage)(ts, stagetime));
3207   PetscFunctionReturn(PETSC_SUCCESS);
3208 }
3209 
3210 /*@
3211   TSPostStage - Runs the user-defined post-stage function set using `TSSetPostStage()`
3212 
3213   Collective
3214 
3215   Input Parameters:
3216 + ts         - The `TS` context obtained from `TSCreate()`
3217 . stagetime  - The absolute time of the current stage
3218 . stageindex - Stage number
3219 - Y          - Array of vectors (of size = total number of stages) with the stage solutions
3220 
3221   Level: developer
3222 
3223   Note:
3224   `TSPostStage()` is typically used within time stepping implementations,
3225   most users would not generally call this routine themselves.
3226 
3227 .seealso: [](ch_ts), `TS`, `TSPreStage()`, `TSSetPreStep()`, `TSPreStep()`, `TSPostStep()`
3228 @*/
3229 PetscErrorCode TSPostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
3230 {
3231   PetscFunctionBegin;
3232   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3233   if (ts->poststage) PetscCallBack("TS callback poststage", (*ts->poststage)(ts, stagetime, stageindex, Y));
3234   PetscFunctionReturn(PETSC_SUCCESS);
3235 }
3236 
3237 /*@
3238   TSPostEvaluate - Runs the user-defined post-evaluate function set using `TSSetPostEvaluate()`
3239 
3240   Collective
3241 
3242   Input Parameter:
3243 . ts - The `TS` context obtained from `TSCreate()`
3244 
3245   Level: developer
3246 
3247   Note:
3248   `TSPostEvaluate()` is typically used within time stepping implementations,
3249   most users would not generally call this routine themselves.
3250 
3251 .seealso: [](ch_ts), `TS`, `TSSetPostEvaluate()`, `TSSetPreStep()`, `TSPreStep()`, `TSPostStep()`
3252 @*/
3253 PetscErrorCode TSPostEvaluate(TS ts)
3254 {
3255   PetscFunctionBegin;
3256   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3257   if (ts->postevaluate) {
3258     Vec              U;
3259     PetscObjectState sprev, spost;
3260 
3261     PetscCall(TSGetSolution(ts, &U));
3262     PetscCall(PetscObjectStateGet((PetscObject)U, &sprev));
3263     PetscCallBack("TS callback postevaluate", (*ts->postevaluate)(ts));
3264     PetscCall(PetscObjectStateGet((PetscObject)U, &spost));
3265     if (sprev != spost) PetscCall(TSRestartStep(ts));
3266   }
3267   PetscFunctionReturn(PETSC_SUCCESS);
3268 }
3269 
3270 /*@C
3271   TSSetPostStep - Sets the general-purpose function
3272   called once at the end of each successful time step.
3273 
3274   Logically Collective
3275 
3276   Input Parameters:
3277 + ts   - The `TS` context obtained from `TSCreate()`
3278 - func - The function
3279 
3280   Calling sequence of `func`:
3281 . ts - the `TS` context
3282 
3283   Level: intermediate
3284 
3285   Note:
3286   The function set by `TSSetPostStep()` is called after each successful step. If the event handler locates an event at the
3287   given step, and `postevent()` modifies the solution vector, the solution vector obtained by `TSGetSolution()` inside `func` will
3288   contain the changes. To get the solution without these changes, use `TSSetPostEvaluate()` to set the appropriate callback.
3289   The scheme of the relevant function calls in `TSSolve()` is shown below
3290 .vb
3291   ...
3292   Step()
3293   PostEvaluate()
3294   EventHandling()
3295   step_rollback ? PostEvaluate() : PostStep()
3296   ...
3297 .ve
3298   where EventHandling() may result in one of the following three outcomes
3299 .vb
3300   (1) | successful step | solution intact
3301   (2) | successful step | solution modified by `postevent()`
3302   (3) | step_rollback   | solution rolled back
3303 .ve
3304 
3305 .seealso: [](ch_ts), `TS`, `TSSetPreStep()`, `TSSetPreStage()`, `TSSetPostEvaluate()`, `TSGetTimeStep()`, `TSGetStepNumber()`, `TSGetTime()`, `TSRestartStep()`
3306 @*/
3307 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS ts))
3308 {
3309   PetscFunctionBegin;
3310   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3311   ts->poststep = func;
3312   PetscFunctionReturn(PETSC_SUCCESS);
3313 }
3314 
3315 /*@
3316   TSPostStep - Runs the user-defined post-step function that was set with `TSSetPostStep()`
3317 
3318   Collective
3319 
3320   Input Parameter:
3321 . ts - The `TS` context obtained from `TSCreate()`
3322 
3323   Note:
3324   `TSPostStep()` is typically used within time stepping implementations,
3325   so most users would not generally call this routine themselves.
3326 
3327   Level: developer
3328 
3329 .seealso: [](ch_ts), `TS`, `TSSetPreStep()`, `TSSetPreStage()`, `TSSetPostEvaluate()`, `TSGetTimeStep()`, `TSGetStepNumber()`, `TSGetTime()`, `TSSetPotsStep()`
3330 @*/
3331 PetscErrorCode TSPostStep(TS ts)
3332 {
3333   PetscFunctionBegin;
3334   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3335   if (ts->poststep) {
3336     Vec              U;
3337     PetscObjectId    idprev;
3338     PetscBool        sameObject;
3339     PetscObjectState sprev, spost;
3340 
3341     PetscCall(TSGetSolution(ts, &U));
3342     PetscCall(PetscObjectGetId((PetscObject)U, &idprev));
3343     PetscCall(PetscObjectStateGet((PetscObject)U, &sprev));
3344     PetscCallBack("TS callback poststep", (*ts->poststep)(ts));
3345     PetscCall(TSGetSolution(ts, &U));
3346     PetscCall(PetscObjectCompareId((PetscObject)U, idprev, &sameObject));
3347     PetscCall(PetscObjectStateGet((PetscObject)U, &spost));
3348     if (!sameObject || sprev != spost) PetscCall(TSRestartStep(ts));
3349   }
3350   PetscFunctionReturn(PETSC_SUCCESS);
3351 }
3352 
3353 /*@
3354   TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
3355 
3356   Collective
3357 
3358   Input Parameters:
3359 + ts - time stepping context
3360 - t  - time to interpolate to
3361 
3362   Output Parameter:
3363 . U - state at given time
3364 
3365   Level: intermediate
3366 
3367   Developer Notes:
3368   `TSInterpolate()` and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
3369 
3370 .seealso: [](ch_ts), `TS`, `TSSetExactFinalTime()`, `TSSolve()`
3371 @*/
3372 PetscErrorCode TSInterpolate(TS ts, PetscReal t, Vec U)
3373 {
3374   PetscFunctionBegin;
3375   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3376   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
3377   PetscCheck(t >= ts->ptime_prev && t <= ts->ptime, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Requested time %g not in last time steps [%g,%g]", (double)t, (double)ts->ptime_prev, (double)ts->ptime);
3378   PetscUseTypeMethod(ts, interpolate, t, U);
3379   PetscFunctionReturn(PETSC_SUCCESS);
3380 }
3381 
3382 /*@
3383   TSStep - Steps one time step
3384 
3385   Collective
3386 
3387   Input Parameter:
3388 . ts - the `TS` context obtained from `TSCreate()`
3389 
3390   Level: developer
3391 
3392   Notes:
3393   The public interface for the ODE/DAE solvers is `TSSolve()`, you should almost for sure be using that routine and not this routine.
3394 
3395   The hook set using `TSSetPreStep()` is called before each attempt to take the step. In general, the time step size may
3396   be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.
3397 
3398   This may over-step the final time provided in `TSSetMaxTime()` depending on the time-step used. `TSSolve()` interpolates to exactly the
3399   time provided in `TSSetMaxTime()`. One can use `TSInterpolate()` to determine an interpolated solution within the final timestep.
3400 
3401 .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetUp()`, `TSDestroy()`, `TSSolve()`, `TSSetPreStep()`, `TSSetPreStage()`, `TSSetPostStage()`, `TSInterpolate()`
3402 @*/
3403 PetscErrorCode TSStep(TS ts)
3404 {
3405   static PetscBool cite = PETSC_FALSE;
3406   PetscReal        ptime;
3407 
3408   PetscFunctionBegin;
3409   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3410   PetscCall(PetscCitationsRegister("@article{tspaper,\n"
3411                                    "  title         = {{PETSc/TS}: A Modern Scalable {DAE/ODE} Solver Library},\n"
3412                                    "  author        = {Abhyankar, Shrirang and Brown, Jed and Constantinescu, Emil and Ghosh, Debojyoti and Smith, Barry F. and Zhang, Hong},\n"
3413                                    "  journal       = {arXiv e-preprints},\n"
3414                                    "  eprint        = {1806.01437},\n"
3415                                    "  archivePrefix = {arXiv},\n"
3416                                    "  year          = {2018}\n}\n",
3417                                    &cite));
3418   PetscCall(TSSetUp(ts));
3419   PetscCall(TSTrajectorySetUp(ts->trajectory, ts));
3420   if (ts->eval_times)
3421     ts->eval_times->worktol = 0; /* In each step of TSSolve() 'eval_times->worktol' will be meaningfully defined (later) only once:
3422                                                    in TSAdaptChoose() or TSEvent_dt_cap(), and then reused till the end of the step */
3423 
3424   PetscCheck(ts->max_time < PETSC_MAX_REAL || ts->max_steps != PETSC_INT_MAX, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "You must call TSSetMaxTime() or TSSetMaxSteps(), or use -ts_max_time <time> or -ts_max_steps <steps>");
3425   PetscCheck(ts->exact_final_time != TS_EXACTFINALTIME_UNSPECIFIED, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "You must call TSSetExactFinalTime() or use -ts_exact_final_time <stepover,interpolate,matchstep> before calling TSStep()");
3426   PetscCheck(ts->exact_final_time != TS_EXACTFINALTIME_MATCHSTEP || ts->adapt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Since TS is not adaptive you cannot use TS_EXACTFINALTIME_MATCHSTEP, suggest TS_EXACTFINALTIME_INTERPOLATE");
3427 
3428   if (!ts->vec_sol0) PetscCall(VecDuplicate(ts->vec_sol, &ts->vec_sol0));
3429   PetscCall(VecCopy(ts->vec_sol, ts->vec_sol0));
3430   ts->time_step0 = ts->time_step;
3431 
3432   if (!ts->steps) ts->ptime_prev = ts->ptime;
3433   ptime = ts->ptime;
3434 
3435   ts->ptime_prev_rollback = ts->ptime_prev;
3436   ts->reason              = TS_CONVERGED_ITERATING;
3437 
3438   PetscCall(PetscLogEventBegin(TS_Step, ts, 0, 0, 0));
3439   PetscUseTypeMethod(ts, step);
3440   PetscCall(PetscLogEventEnd(TS_Step, ts, 0, 0, 0));
3441 
3442   if (ts->reason >= 0) {
3443     ts->ptime_prev = ptime;
3444     ts->steps++;
3445     ts->steprollback = PETSC_FALSE;
3446     ts->steprestart  = PETSC_FALSE;
3447     ts->stepresize   = PETSC_FALSE;
3448   }
3449 
3450   if (ts->reason < 0 && ts->errorifstepfailed) {
3451     PetscCall(TSMonitorCancel(ts));
3452     PetscCheck(ts->reason != TS_DIVERGED_NONLINEAR_SOLVE, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSStep has failed due to %s, increase -ts_max_snes_failures or use unlimited to attempt recovery", TSConvergedReasons[ts->reason]);
3453     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSStep has failed due to %s", TSConvergedReasons[ts->reason]);
3454   }
3455   PetscFunctionReturn(PETSC_SUCCESS);
3456 }
3457 
3458 /*@
3459   TSEvaluateWLTE - Evaluate the weighted local truncation error norm
3460   at the end of a time step with a given order of accuracy.
3461 
3462   Collective
3463 
3464   Input Parameters:
3465 + ts        - time stepping context
3466 - wnormtype - norm type, either `NORM_2` or `NORM_INFINITY`
3467 
3468   Input/Output Parameter:
3469 . order - optional, desired order for the error evaluation or `PETSC_DECIDE`;
3470            on output, the actual order of the error evaluation
3471 
3472   Output Parameter:
3473 . wlte - the weighted local truncation error norm
3474 
3475   Level: advanced
3476 
3477   Note:
3478   If the timestepper cannot evaluate the error in a particular step
3479   (eg. in the first step or restart steps after event handling),
3480   this routine returns wlte=-1.0 .
3481 
3482 .seealso: [](ch_ts), `TS`, `TSStep()`, `TSAdapt`, `TSErrorWeightedNorm()`
3483 @*/
3484 PetscErrorCode TSEvaluateWLTE(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
3485 {
3486   PetscFunctionBegin;
3487   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3488   PetscValidType(ts, 1);
3489   PetscValidLogicalCollectiveEnum(ts, wnormtype, 2);
3490   if (order) PetscAssertPointer(order, 3);
3491   if (order) PetscValidLogicalCollectiveInt(ts, *order, 3);
3492   PetscAssertPointer(wlte, 4);
3493   PetscCheck(wnormtype == NORM_2 || wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "No support for norm type %s", NormTypes[wnormtype]);
3494   PetscUseTypeMethod(ts, evaluatewlte, wnormtype, order, wlte);
3495   PetscFunctionReturn(PETSC_SUCCESS);
3496 }
3497 
3498 /*@
3499   TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
3500 
3501   Collective
3502 
3503   Input Parameters:
3504 + ts    - time stepping context
3505 . order - desired order of accuracy
3506 - done  - whether the step was evaluated at this order (pass `NULL` to generate an error if not available)
3507 
3508   Output Parameter:
3509 . U - state at the end of the current step
3510 
3511   Level: advanced
3512 
3513   Notes:
3514   This function cannot be called until all stages have been evaluated.
3515 
3516   It is normally called by adaptive controllers before a step has been accepted and may also be called by the user after `TSStep()` has returned.
3517 
3518 .seealso: [](ch_ts), `TS`, `TSStep()`, `TSAdapt`
3519 @*/
3520 PetscErrorCode TSEvaluateStep(TS ts, PetscInt order, Vec U, PetscBool *done)
3521 {
3522   PetscFunctionBegin;
3523   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3524   PetscValidType(ts, 1);
3525   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
3526   PetscUseTypeMethod(ts, evaluatestep, order, U, done);
3527   PetscFunctionReturn(PETSC_SUCCESS);
3528 }
3529 
3530 /*@C
3531   TSGetComputeInitialCondition - Get the function used to automatically compute an initial condition for the timestepping.
3532 
3533   Not collective
3534 
3535   Input Parameter:
3536 . ts - time stepping context
3537 
3538   Output Parameter:
3539 . initCondition - The function which computes an initial condition
3540 
3541   Calling sequence of `initCondition`:
3542 + ts - The timestepping context
3543 - u  - The input vector in which the initial condition is stored
3544 
3545   Level: advanced
3546 
3547 .seealso: [](ch_ts), `TS`, `TSSetComputeInitialCondition()`, `TSComputeInitialCondition()`
3548 @*/
3549 PetscErrorCode TSGetComputeInitialCondition(TS ts, PetscErrorCode (**initCondition)(TS ts, Vec u))
3550 {
3551   PetscFunctionBegin;
3552   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3553   PetscAssertPointer(initCondition, 2);
3554   *initCondition = ts->ops->initcondition;
3555   PetscFunctionReturn(PETSC_SUCCESS);
3556 }
3557 
3558 /*@C
3559   TSSetComputeInitialCondition - Set the function used to automatically compute an initial condition for the timestepping.
3560 
3561   Logically collective
3562 
3563   Input Parameters:
3564 + ts            - time stepping context
3565 - initCondition - The function which computes an initial condition
3566 
3567   Calling sequence of `initCondition`:
3568 + ts - The timestepping context
3569 - e  - The input vector in which the initial condition is to be stored
3570 
3571   Level: advanced
3572 
3573 .seealso: [](ch_ts), `TS`, `TSGetComputeInitialCondition()`, `TSComputeInitialCondition()`
3574 @*/
3575 PetscErrorCode TSSetComputeInitialCondition(TS ts, PetscErrorCode (*initCondition)(TS ts, Vec e))
3576 {
3577   PetscFunctionBegin;
3578   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3579   PetscValidFunction(initCondition, 2);
3580   ts->ops->initcondition = initCondition;
3581   PetscFunctionReturn(PETSC_SUCCESS);
3582 }
3583 
3584 /*@
3585   TSComputeInitialCondition - Compute an initial condition for the timestepping using the function previously set with `TSSetComputeInitialCondition()`
3586 
3587   Collective
3588 
3589   Input Parameters:
3590 + ts - time stepping context
3591 - u  - The `Vec` to store the condition in which will be used in `TSSolve()`
3592 
3593   Level: advanced
3594 
3595 .seealso: [](ch_ts), `TS`, `TSGetComputeInitialCondition()`, `TSSetComputeInitialCondition()`, `TSSolve()`
3596 @*/
3597 PetscErrorCode TSComputeInitialCondition(TS ts, Vec u)
3598 {
3599   PetscFunctionBegin;
3600   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3601   PetscValidHeaderSpecific(u, VEC_CLASSID, 2);
3602   PetscTryTypeMethod(ts, initcondition, u);
3603   PetscFunctionReturn(PETSC_SUCCESS);
3604 }
3605 
3606 /*@C
3607   TSGetComputeExactError - Get the function used to automatically compute the exact error for the timestepping.
3608 
3609   Not collective
3610 
3611   Input Parameter:
3612 . ts - time stepping context
3613 
3614   Output Parameter:
3615 . exactError - The function which computes the solution error
3616 
3617   Calling sequence of `exactError`:
3618 + ts - The timestepping context
3619 . u  - The approximate solution vector
3620 - e  - The vector in which the error is stored
3621 
3622   Level: advanced
3623 
3624 .seealso: [](ch_ts), `TS`, `TSComputeExactError()`
3625 @*/
3626 PetscErrorCode TSGetComputeExactError(TS ts, PetscErrorCode (**exactError)(TS ts, Vec u, Vec e))
3627 {
3628   PetscFunctionBegin;
3629   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3630   PetscAssertPointer(exactError, 2);
3631   *exactError = ts->ops->exacterror;
3632   PetscFunctionReturn(PETSC_SUCCESS);
3633 }
3634 
3635 /*@C
3636   TSSetComputeExactError - Set the function used to automatically compute the exact error for the timestepping.
3637 
3638   Logically collective
3639 
3640   Input Parameters:
3641 + ts         - time stepping context
3642 - exactError - The function which computes the solution error
3643 
3644   Calling sequence of `exactError`:
3645 + ts - The timestepping context
3646 . u  - The approximate solution vector
3647 - e  - The  vector in which the error is stored
3648 
3649   Level: advanced
3650 
3651 .seealso: [](ch_ts), `TS`, `TSGetComputeExactError()`, `TSComputeExactError()`
3652 @*/
3653 PetscErrorCode TSSetComputeExactError(TS ts, PetscErrorCode (*exactError)(TS ts, Vec u, Vec e))
3654 {
3655   PetscFunctionBegin;
3656   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3657   PetscValidFunction(exactError, 2);
3658   ts->ops->exacterror = exactError;
3659   PetscFunctionReturn(PETSC_SUCCESS);
3660 }
3661 
3662 /*@
3663   TSComputeExactError - Compute the solution error for the timestepping using the function previously set with `TSSetComputeExactError()`
3664 
3665   Collective
3666 
3667   Input Parameters:
3668 + ts - time stepping context
3669 . u  - The approximate solution
3670 - e  - The `Vec` used to store the error
3671 
3672   Level: advanced
3673 
3674 .seealso: [](ch_ts), `TS`, `TSGetComputeInitialCondition()`, `TSSetComputeInitialCondition()`, `TSSolve()`
3675 @*/
3676 PetscErrorCode TSComputeExactError(TS ts, Vec u, Vec e)
3677 {
3678   PetscFunctionBegin;
3679   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3680   PetscValidHeaderSpecific(u, VEC_CLASSID, 2);
3681   PetscValidHeaderSpecific(e, VEC_CLASSID, 3);
3682   PetscTryTypeMethod(ts, exacterror, u, e);
3683   PetscFunctionReturn(PETSC_SUCCESS);
3684 }
3685 
3686 /*@C
3687   TSSetResize - Sets the resize callbacks.
3688 
3689   Logically Collective
3690 
3691   Input Parameters:
3692 + ts       - The `TS` context obtained from `TSCreate()`
3693 . rollback - Whether a resize will restart the step
3694 . setup    - The setup function
3695 . transfer - The transfer function
3696 - ctx      - [optional] The user-defined context
3697 
3698   Calling sequence of `setup`:
3699 + ts     - the `TS` context
3700 . step   - the current step
3701 . time   - the current time
3702 . state  - the current vector of state
3703 . resize - (output parameter) `PETSC_TRUE` if need resizing, `PETSC_FALSE` otherwise
3704 - ctx    - user defined context
3705 
3706   Calling sequence of `transfer`:
3707 + ts      - the `TS` context
3708 . nv      - the number of vectors to be transferred
3709 . vecsin  - array of vectors to be transferred
3710 . vecsout - array of transferred vectors
3711 - ctx     - user defined context
3712 
3713   Notes:
3714   The `setup` function is called inside `TSSolve()` after `TSEventHandler()` or after `TSPostStep()`
3715   depending on the `rollback` value: if `rollback` is true, then these callbacks behave as error indicators
3716   and will flag the need to remesh and restart the current step. Otherwise, they will simply flag the solver
3717   that the size of the discrete problem has changed.
3718   In both cases, the solver will collect the needed vectors that will be
3719   transferred from the old to the new sizes using the `transfer` callback. These vectors will include the
3720   current solution vector, and other vectors needed by the specific solver used.
3721   For example, `TSBDF` uses previous solutions vectors to solve for the next time step.
3722   Other application specific objects associated with the solver, i.e. Jacobian matrices and `DM`,
3723   will be automatically reset if the sizes are changed and they must be specified again by the user
3724   inside the `transfer` function.
3725   The input and output arrays passed to `transfer` are allocated by PETSc.
3726   Vectors in `vecsout` must be created by the user.
3727   Ownership of vectors in `vecsout` is transferred to PETSc.
3728 
3729   Level: advanced
3730 
3731 .seealso: [](ch_ts), `TS`, `TSSetDM()`, `TSSetIJacobian()`, `TSSetRHSJacobian()`
3732 @*/
3733 PetscErrorCode TSSetResize(TS ts, PetscBool rollback, PetscErrorCode (*setup)(TS ts, PetscInt step, PetscReal time, Vec state, PetscBool *resize, void *ctx), PetscErrorCode (*transfer)(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], void *ctx), void *ctx)
3734 {
3735   PetscFunctionBegin;
3736   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3737   ts->resizerollback = rollback;
3738   ts->resizesetup    = setup;
3739   ts->resizetransfer = transfer;
3740   ts->resizectx      = ctx;
3741   PetscFunctionReturn(PETSC_SUCCESS);
3742 }
3743 
3744 /*
3745   TSResizeRegisterOrRetrieve - Register or import vectors transferred with `TSResize()`.
3746 
3747   Collective
3748 
3749   Input Parameters:
3750 + ts   - The `TS` context obtained from `TSCreate()`
3751 - flg - If `PETSC_TRUE` each TS implementation (e.g. `TSBDF`) will register vectors to be transferred, if `PETSC_FALSE` vectors will be imported from transferred vectors.
3752 
3753   Level: developer
3754 
3755   Note:
3756   `TSResizeRegisterOrRetrieve()` is declared PETSC_INTERN since it is
3757    used within time stepping implementations,
3758    so most users would not generally call this routine themselves.
3759 
3760 .seealso: [](ch_ts), `TS`, `TSSetResize()`
3761 @*/
3762 static PetscErrorCode TSResizeRegisterOrRetrieve(TS ts, PetscBool flg)
3763 {
3764   PetscFunctionBegin;
3765   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3766   PetscTryTypeMethod(ts, resizeregister, flg);
3767   /* PetscTryTypeMethod(adapt, resizeregister, flg); */
3768   PetscFunctionReturn(PETSC_SUCCESS);
3769 }
3770 
3771 static PetscErrorCode TSResizeReset(TS ts)
3772 {
3773   PetscFunctionBegin;
3774   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3775   PetscCall(PetscObjectListDestroy(&ts->resizetransferobjs));
3776   PetscFunctionReturn(PETSC_SUCCESS);
3777 }
3778 
3779 static PetscErrorCode TSResizeTransferVecs(TS ts, PetscInt cnt, Vec vecsin[], Vec vecsout[])
3780 {
3781   PetscFunctionBegin;
3782   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3783   PetscValidLogicalCollectiveInt(ts, cnt, 2);
3784   for (PetscInt i = 0; i < cnt; i++) PetscCall(VecLockReadPush(vecsin[i]));
3785   if (ts->resizetransfer) {
3786     PetscCall(PetscInfo(ts, "Transferring %" PetscInt_FMT " vectors\n", cnt));
3787     PetscCallBack("TS callback resize transfer", (*ts->resizetransfer)(ts, cnt, vecsin, vecsout, ts->resizectx));
3788   }
3789   for (PetscInt i = 0; i < cnt; i++) PetscCall(VecLockReadPop(vecsin[i]));
3790   PetscFunctionReturn(PETSC_SUCCESS);
3791 }
3792 
3793 /*@C
3794   TSResizeRegisterVec - Register a vector to be transferred with `TSResize()`.
3795 
3796   Collective
3797 
3798   Input Parameters:
3799 + ts   - The `TS` context obtained from `TSCreate()`
3800 . name - A string identifying the vector
3801 - vec  - The vector
3802 
3803   Level: developer
3804 
3805   Note:
3806   `TSResizeRegisterVec()` is typically used within time stepping implementations,
3807   so most users would not generally call this routine themselves.
3808 
3809 .seealso: [](ch_ts), `TS`, `TSSetResize()`, `TSResize()`, `TSResizeRetrieveVec()`
3810 @*/
3811 PetscErrorCode TSResizeRegisterVec(TS ts, const char name[], Vec vec)
3812 {
3813   PetscFunctionBegin;
3814   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3815   PetscAssertPointer(name, 2);
3816   if (vec) PetscValidHeaderSpecific(vec, VEC_CLASSID, 3);
3817   PetscCall(PetscObjectListAdd(&ts->resizetransferobjs, name, (PetscObject)vec));
3818   PetscFunctionReturn(PETSC_SUCCESS);
3819 }
3820 
3821 /*@C
3822   TSResizeRetrieveVec - Retrieve a vector registered with `TSResizeRegisterVec()`.
3823 
3824   Collective
3825 
3826   Input Parameters:
3827 + ts   - The `TS` context obtained from `TSCreate()`
3828 . name - A string identifying the vector
3829 - vec  - The vector
3830 
3831   Level: developer
3832 
3833   Note:
3834   `TSResizeRetrieveVec()` is typically used within time stepping implementations,
3835   so most users would not generally call this routine themselves.
3836 
3837 .seealso: [](ch_ts), `TS`, `TSSetResize()`, `TSResize()`, `TSResizeRegisterVec()`
3838 @*/
3839 PetscErrorCode TSResizeRetrieveVec(TS ts, const char name[], Vec *vec)
3840 {
3841   PetscFunctionBegin;
3842   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3843   PetscAssertPointer(name, 2);
3844   PetscAssertPointer(vec, 3);
3845   PetscCall(PetscObjectListFind(ts->resizetransferobjs, name, (PetscObject *)vec));
3846   PetscFunctionReturn(PETSC_SUCCESS);
3847 }
3848 
3849 static PetscErrorCode TSResizeGetVecArray(TS ts, PetscInt *nv, const char **names[], Vec *vecs[])
3850 {
3851   PetscInt        cnt;
3852   PetscObjectList tmp;
3853   Vec            *vecsin  = NULL;
3854   const char    **namesin = NULL;
3855 
3856   PetscFunctionBegin;
3857   for (tmp = ts->resizetransferobjs, cnt = 0; tmp; tmp = tmp->next)
3858     if (tmp->obj && tmp->obj->classid == VEC_CLASSID) cnt++;
3859   if (names) PetscCall(PetscMalloc1(cnt, &namesin));
3860   if (vecs) PetscCall(PetscMalloc1(cnt, &vecsin));
3861   for (tmp = ts->resizetransferobjs, cnt = 0; tmp; tmp = tmp->next) {
3862     if (tmp->obj && tmp->obj->classid == VEC_CLASSID) {
3863       if (vecs) vecsin[cnt] = (Vec)tmp->obj;
3864       if (names) namesin[cnt] = tmp->name;
3865       cnt++;
3866     }
3867   }
3868   if (nv) *nv = cnt;
3869   if (names) *names = namesin;
3870   if (vecs) *vecs = vecsin;
3871   PetscFunctionReturn(PETSC_SUCCESS);
3872 }
3873 
3874 /*@
3875   TSResize - Runs the user-defined transfer functions provided with `TSSetResize()`
3876 
3877   Collective
3878 
3879   Input Parameter:
3880 . ts - The `TS` context obtained from `TSCreate()`
3881 
3882   Level: developer
3883 
3884   Note:
3885   `TSResize()` is typically used within time stepping implementations,
3886   so most users would not generally call this routine themselves.
3887 
3888 .seealso: [](ch_ts), `TS`, `TSSetResize()`
3889 @*/
3890 PetscErrorCode TSResize(TS ts)
3891 {
3892   PetscInt     nv      = 0;
3893   const char **names   = NULL;
3894   Vec         *vecsin  = NULL;
3895   const char  *solname = "ts:vec_sol";
3896 
3897   PetscFunctionBegin;
3898   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3899   if (!ts->resizesetup) PetscFunctionReturn(PETSC_SUCCESS);
3900   if (ts->resizesetup) {
3901     PetscCall(VecLockReadPush(ts->vec_sol));
3902     PetscCallBack("TS callback resize setup", (*ts->resizesetup)(ts, ts->steps, ts->ptime, ts->vec_sol, &ts->stepresize, ts->resizectx));
3903     PetscCall(VecLockReadPop(ts->vec_sol));
3904     if (ts->stepresize) {
3905       if (ts->resizerollback) {
3906         PetscCall(TSRollBack(ts));
3907         ts->time_step = ts->time_step0;
3908       }
3909       PetscCall(TSResizeRegisterVec(ts, solname, ts->vec_sol));
3910       PetscCall(TSResizeRegisterOrRetrieve(ts, PETSC_TRUE)); /* specific impls register their own objects */
3911     }
3912   }
3913 
3914   PetscCall(TSResizeGetVecArray(ts, &nv, &names, &vecsin));
3915   if (nv) {
3916     Vec *vecsout, vecsol;
3917 
3918     /* Reset internal objects */
3919     PetscCall(TSReset(ts));
3920 
3921     /* Transfer needed vectors (users can call SetJacobian, SetDM, etc. here) */
3922     PetscCall(PetscCalloc1(nv, &vecsout));
3923     PetscCall(TSResizeTransferVecs(ts, nv, vecsin, vecsout));
3924     for (PetscInt i = 0; i < nv; i++) {
3925       const char *name;
3926       char       *oname;
3927 
3928       PetscCall(PetscObjectGetName((PetscObject)vecsin[i], &name));
3929       PetscCall(PetscStrallocpy(name, &oname));
3930       PetscCall(TSResizeRegisterVec(ts, names[i], vecsout[i]));
3931       if (vecsout[i]) PetscCall(PetscObjectSetName((PetscObject)vecsout[i], oname));
3932       PetscCall(PetscFree(oname));
3933       PetscCall(VecDestroy(&vecsout[i]));
3934     }
3935     PetscCall(PetscFree(vecsout));
3936     PetscCall(TSResizeRegisterOrRetrieve(ts, PETSC_FALSE)); /* specific impls import the transferred objects */
3937 
3938     PetscCall(TSResizeRetrieveVec(ts, solname, &vecsol));
3939     if (vecsol) PetscCall(TSSetSolution(ts, vecsol));
3940     PetscAssert(ts->vec_sol, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_NULL, "Missing TS solution");
3941   }
3942 
3943   PetscCall(PetscFree(names));
3944   PetscCall(PetscFree(vecsin));
3945   PetscCall(TSResizeReset(ts));
3946   PetscFunctionReturn(PETSC_SUCCESS);
3947 }
3948 
3949 /*@
3950   TSSolve - Steps the requested number of timesteps.
3951 
3952   Collective
3953 
3954   Input Parameters:
3955 + ts - the `TS` context obtained from `TSCreate()`
3956 - u  - the solution vector  (can be null if `TSSetSolution()` was used and `TSSetExactFinalTime`(ts,`TS_EXACTFINALTIME_MATCHSTEP`) was not used,
3957        otherwise it must contain the initial conditions and will contain the solution at the final requested time
3958 
3959   Level: beginner
3960 
3961   Notes:
3962   The final time returned by this function may be different from the time of the internally
3963   held state accessible by `TSGetSolution()` and `TSGetTime()` because the method may have
3964   stepped over the final time.
3965 
3966 .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetSolution()`, `TSStep()`, `TSGetTime()`, `TSGetSolveTime()`
3967 @*/
3968 PetscErrorCode TSSolve(TS ts, Vec u)
3969 {
3970   Vec solution;
3971 
3972   PetscFunctionBegin;
3973   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3974   if (u) PetscValidHeaderSpecific(u, VEC_CLASSID, 2);
3975 
3976   PetscCall(TSSetExactFinalTimeDefault(ts));
3977   if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && u) { /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
3978     if (!ts->vec_sol || u == ts->vec_sol) {
3979       PetscCall(VecDuplicate(u, &solution));
3980       PetscCall(TSSetSolution(ts, solution));
3981       PetscCall(VecDestroy(&solution)); /* grant ownership */
3982     }
3983     PetscCall(VecCopy(u, ts->vec_sol));
3984     PetscCheck(!ts->forward_solve, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Sensitivity analysis does not support the mode TS_EXACTFINALTIME_INTERPOLATE");
3985   } else if (u) PetscCall(TSSetSolution(ts, u));
3986   PetscCall(TSSetUp(ts));
3987   PetscCall(TSTrajectorySetUp(ts->trajectory, ts));
3988 
3989   PetscCheck(ts->max_time < PETSC_MAX_REAL || ts->max_steps != PETSC_INT_MAX, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "You must call TSSetMaxTime() or TSSetMaxSteps(), or use -ts_max_time <time> or -ts_max_steps <steps>");
3990   PetscCheck(ts->exact_final_time != TS_EXACTFINALTIME_UNSPECIFIED, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "You must call TSSetExactFinalTime() or use -ts_exact_final_time <stepover,interpolate,matchstep> before calling TSSolve()");
3991   PetscCheck(ts->exact_final_time != TS_EXACTFINALTIME_MATCHSTEP || ts->adapt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Since TS is not adaptive you cannot use TS_EXACTFINALTIME_MATCHSTEP, suggest TS_EXACTFINALTIME_INTERPOLATE");
3992   PetscCheck(!(ts->eval_times && ts->exact_final_time != TS_EXACTFINALTIME_MATCHSTEP), PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "You must use TS_EXACTFINALTIME_MATCHSTEP when using time span or evaluation times");
3993 
3994   if (ts->eval_times) {
3995     for (PetscInt i = 0; i < ts->eval_times->num_time_points; i++) {
3996       PetscBool is_close = PetscIsCloseAtTol(ts->ptime, ts->eval_times->time_points[i], ts->eval_times->reltol * ts->time_step + ts->eval_times->abstol, 0);
3997       if (ts->ptime <= ts->eval_times->time_points[i] || is_close) {
3998         ts->eval_times->time_point_idx = i;
3999         if (is_close) { /* starting point in evaluation times */
4000           PetscCall(VecCopy(ts->vec_sol, ts->eval_times->sol_vecs[ts->eval_times->sol_ctr]));
4001           ts->eval_times->sol_times[ts->eval_times->sol_ctr] = ts->ptime;
4002           ts->eval_times->sol_ctr++;
4003           ts->eval_times->time_point_idx++;
4004         }
4005         break;
4006       }
4007     }
4008   }
4009 
4010   if (ts->forward_solve) PetscCall(TSForwardSetUp(ts));
4011 
4012   /* reset number of steps only when the step is not restarted. ARKIMEX
4013      restarts the step after an event. Resetting these counters in such case causes
4014      TSTrajectory to incorrectly save the output files
4015   */
4016   /* reset time step and iteration counters */
4017   if (!ts->steps) {
4018     ts->ksp_its           = 0;
4019     ts->snes_its          = 0;
4020     ts->num_snes_failures = 0;
4021     ts->reject            = 0;
4022     ts->steprestart       = PETSC_TRUE;
4023     ts->steprollback      = PETSC_FALSE;
4024     ts->stepresize        = PETSC_FALSE;
4025     ts->rhsjacobian.time  = PETSC_MIN_REAL;
4026   }
4027 
4028   /* make sure initial time step does not overshoot final time or the next point in evaluation times */
4029   if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
4030     PetscReal maxdt;
4031     PetscReal dt = ts->time_step;
4032 
4033     if (ts->eval_times) maxdt = ts->eval_times->time_points[ts->eval_times->time_point_idx] - ts->ptime;
4034     else maxdt = ts->max_time - ts->ptime;
4035     ts->time_step = dt >= maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt);
4036   }
4037   ts->reason = TS_CONVERGED_ITERATING;
4038 
4039   {
4040     PetscViewer       viewer;
4041     PetscViewerFormat format;
4042     PetscBool         flg;
4043     static PetscBool  incall = PETSC_FALSE;
4044 
4045     if (!incall) {
4046       /* Estimate the convergence rate of the time discretization */
4047       PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_convergence_estimate", &viewer, &format, &flg));
4048       if (flg) {
4049         PetscConvEst conv;
4050         DM           dm;
4051         PetscReal   *alpha; /* Convergence rate of the solution error for each field in the L_2 norm */
4052         PetscInt     Nf;
4053         PetscBool    checkTemporal = PETSC_TRUE;
4054 
4055         incall = PETSC_TRUE;
4056         PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_convergence_temporal", &checkTemporal, &flg));
4057         PetscCall(TSGetDM(ts, &dm));
4058         PetscCall(DMGetNumFields(dm, &Nf));
4059         PetscCall(PetscCalloc1(PetscMax(Nf, 1), &alpha));
4060         PetscCall(PetscConvEstCreate(PetscObjectComm((PetscObject)ts), &conv));
4061         PetscCall(PetscConvEstUseTS(conv, checkTemporal));
4062         PetscCall(PetscConvEstSetSolver(conv, (PetscObject)ts));
4063         PetscCall(PetscConvEstSetFromOptions(conv));
4064         PetscCall(PetscConvEstSetUp(conv));
4065         PetscCall(PetscConvEstGetConvRate(conv, alpha));
4066         PetscCall(PetscViewerPushFormat(viewer, format));
4067         PetscCall(PetscConvEstRateView(conv, alpha, viewer));
4068         PetscCall(PetscViewerPopFormat(viewer));
4069         PetscCall(PetscViewerDestroy(&viewer));
4070         PetscCall(PetscConvEstDestroy(&conv));
4071         PetscCall(PetscFree(alpha));
4072         incall = PETSC_FALSE;
4073       }
4074     }
4075   }
4076 
4077   PetscCall(TSViewFromOptions(ts, NULL, "-ts_view_pre"));
4078 
4079   if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */
4080     PetscUseTypeMethod(ts, solve);
4081     if (u) PetscCall(VecCopy(ts->vec_sol, u));
4082     ts->solvetime = ts->ptime;
4083     solution      = ts->vec_sol;
4084   } else { /* Step the requested number of timesteps. */
4085     if (ts->steps >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
4086     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
4087 
4088     if (!ts->steps) {
4089       PetscCall(TSTrajectorySet(ts->trajectory, ts, ts->steps, ts->ptime, ts->vec_sol));
4090       PetscCall(TSEventInitialize(ts->event, ts, ts->ptime, ts->vec_sol));
4091     }
4092 
4093     while (!ts->reason) {
4094       PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol));
4095       if (!ts->steprollback || (ts->stepresize && ts->resizerollback)) PetscCall(TSPreStep(ts));
4096       PetscCall(TSStep(ts));
4097       if (ts->testjacobian) PetscCall(TSRHSJacobianTest(ts, NULL));
4098       if (ts->testjacobiantranspose) PetscCall(TSRHSJacobianTestTranspose(ts, NULL));
4099       if (ts->quadraturets && ts->costintegralfwd) { /* Must evaluate the cost integral before event is handled. The cost integral value can also be rolled back. */
4100         if (ts->reason >= 0) ts->steps--;            /* Revert the step number changed by TSStep() */
4101         PetscCall(TSForwardCostIntegral(ts));
4102         if (ts->reason >= 0) ts->steps++;
4103       }
4104       if (ts->forward_solve) {            /* compute forward sensitivities before event handling because postevent() may change RHS and jump conditions may have to be applied */
4105         if (ts->reason >= 0) ts->steps--; /* Revert the step number changed by TSStep() */
4106         PetscCall(TSForwardStep(ts));
4107         if (ts->reason >= 0) ts->steps++;
4108       }
4109       PetscCall(TSPostEvaluate(ts));
4110       PetscCall(TSEventHandler(ts)); /* The right-hand side may be changed due to event. Be careful with Any computation using the RHS information after this point. */
4111       if (ts->steprollback) PetscCall(TSPostEvaluate(ts));
4112       if (!ts->steprollback && ts->resizerollback) PetscCall(TSResize(ts));
4113       /* check convergence */
4114       if (!ts->reason) {
4115         if (ts->steps >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
4116         else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
4117       }
4118       if (!ts->steprollback) {
4119         PetscCall(TSTrajectorySet(ts->trajectory, ts, ts->steps, ts->ptime, ts->vec_sol));
4120         PetscCall(TSPostStep(ts));
4121         if (!ts->resizerollback) PetscCall(TSResize(ts));
4122 
4123         if (ts->eval_times && ts->eval_times->time_point_idx < ts->eval_times->num_time_points) {
4124           PetscCheck(ts->eval_times->worktol > 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Unexpected state !(eval_times->worktol > 0) in TSSolve()");
4125           if (PetscIsCloseAtTol(ts->ptime, ts->eval_times->time_points[ts->eval_times->time_point_idx], ts->eval_times->worktol, 0)) {
4126             ts->eval_times->sol_times[ts->eval_times->sol_ctr] = ts->ptime;
4127             PetscCall(VecCopy(ts->vec_sol, ts->eval_times->sol_vecs[ts->eval_times->sol_ctr]));
4128             ts->eval_times->sol_ctr++;
4129             ts->eval_times->time_point_idx++;
4130           }
4131         }
4132       }
4133     }
4134     PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol));
4135 
4136     if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && ts->ptime > ts->max_time) {
4137       if (!u) u = ts->vec_sol;
4138       PetscCall(TSInterpolate(ts, ts->max_time, u));
4139       ts->solvetime = ts->max_time;
4140       solution      = u;
4141       PetscCall(TSMonitor(ts, -1, ts->solvetime, solution));
4142     } else {
4143       if (u) PetscCall(VecCopy(ts->vec_sol, u));
4144       ts->solvetime = ts->ptime;
4145       solution      = ts->vec_sol;
4146     }
4147   }
4148 
4149   PetscCall(TSViewFromOptions(ts, NULL, "-ts_view"));
4150   PetscCall(VecViewFromOptions(solution, (PetscObject)ts, "-ts_view_solution"));
4151   PetscCall(PetscObjectSAWsBlock((PetscObject)ts));
4152   if (ts->adjoint_solve) PetscCall(TSAdjointSolve(ts));
4153   PetscFunctionReturn(PETSC_SUCCESS);
4154 }
4155 
4156 /*@
4157   TSGetTime - Gets the time of the most recently completed step.
4158 
4159   Not Collective
4160 
4161   Input Parameter:
4162 . ts - the `TS` context obtained from `TSCreate()`
4163 
4164   Output Parameter:
4165 . t - the current time. This time may not corresponds to the final time set with `TSSetMaxTime()`, use `TSGetSolveTime()`.
4166 
4167   Level: beginner
4168 
4169   Note:
4170   When called during time step evaluation (e.g. during residual evaluation or via hooks set using `TSSetPreStep()`,
4171   `TSSetPreStage()`, `TSSetPostStage()`, or `TSSetPostStep()`), the time is the time at the start of the step being evaluated.
4172 
4173 .seealso: [](ch_ts), `TS`, `TSGetSolveTime()`, `TSSetTime()`, `TSGetTimeStep()`, `TSGetStepNumber()`
4174 @*/
4175 PetscErrorCode TSGetTime(TS ts, PetscReal *t)
4176 {
4177   PetscFunctionBegin;
4178   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4179   PetscAssertPointer(t, 2);
4180   *t = ts->ptime;
4181   PetscFunctionReturn(PETSC_SUCCESS);
4182 }
4183 
4184 /*@
4185   TSGetPrevTime - Gets the starting time of the previously completed step.
4186 
4187   Not Collective
4188 
4189   Input Parameter:
4190 . ts - the `TS` context obtained from `TSCreate()`
4191 
4192   Output Parameter:
4193 . t - the previous time
4194 
4195   Level: beginner
4196 
4197 .seealso: [](ch_ts), `TS`, `TSGetTime()`, `TSGetSolveTime()`, `TSGetTimeStep()`
4198 @*/
4199 PetscErrorCode TSGetPrevTime(TS ts, PetscReal *t)
4200 {
4201   PetscFunctionBegin;
4202   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4203   PetscAssertPointer(t, 2);
4204   *t = ts->ptime_prev;
4205   PetscFunctionReturn(PETSC_SUCCESS);
4206 }
4207 
4208 /*@
4209   TSSetTime - Allows one to reset the time.
4210 
4211   Logically Collective
4212 
4213   Input Parameters:
4214 + ts - the `TS` context obtained from `TSCreate()`
4215 - t  - the time
4216 
4217   Level: intermediate
4218 
4219 .seealso: [](ch_ts), `TS`, `TSGetTime()`, `TSSetMaxSteps()`
4220 @*/
4221 PetscErrorCode TSSetTime(TS ts, PetscReal t)
4222 {
4223   PetscFunctionBegin;
4224   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4225   PetscValidLogicalCollectiveReal(ts, t, 2);
4226   ts->ptime = t;
4227   PetscFunctionReturn(PETSC_SUCCESS);
4228 }
4229 
4230 /*@
4231   TSSetOptionsPrefix - Sets the prefix used for searching for all
4232   TS options in the database.
4233 
4234   Logically Collective
4235 
4236   Input Parameters:
4237 + ts     - The `TS` context
4238 - prefix - The prefix to prepend to all option names
4239 
4240   Level: advanced
4241 
4242   Note:
4243   A hyphen (-) must NOT be given at the beginning of the prefix name.
4244   The first character of all runtime options is AUTOMATICALLY the
4245   hyphen.
4246 
4247 .seealso: [](ch_ts), `TS`, `TSSetFromOptions()`, `TSAppendOptionsPrefix()`
4248 @*/
4249 PetscErrorCode TSSetOptionsPrefix(TS ts, const char prefix[])
4250 {
4251   SNES snes;
4252 
4253   PetscFunctionBegin;
4254   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4255   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ts, prefix));
4256   PetscCall(TSGetSNES(ts, &snes));
4257   PetscCall(SNESSetOptionsPrefix(snes, prefix));
4258   PetscFunctionReturn(PETSC_SUCCESS);
4259 }
4260 
4261 /*@
4262   TSAppendOptionsPrefix - Appends to the prefix used for searching for all
4263   TS options in the database.
4264 
4265   Logically Collective
4266 
4267   Input Parameters:
4268 + ts     - The `TS` context
4269 - prefix - The prefix to prepend to all option names
4270 
4271   Level: advanced
4272 
4273   Note:
4274   A hyphen (-) must NOT be given at the beginning of the prefix name.
4275   The first character of all runtime options is AUTOMATICALLY the
4276   hyphen.
4277 
4278 .seealso: [](ch_ts), `TS`, `TSGetOptionsPrefix()`, `TSSetOptionsPrefix()`, `TSSetFromOptions()`
4279 @*/
4280 PetscErrorCode TSAppendOptionsPrefix(TS ts, const char prefix[])
4281 {
4282   SNES snes;
4283 
4284   PetscFunctionBegin;
4285   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4286   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)ts, prefix));
4287   PetscCall(TSGetSNES(ts, &snes));
4288   PetscCall(SNESAppendOptionsPrefix(snes, prefix));
4289   PetscFunctionReturn(PETSC_SUCCESS);
4290 }
4291 
4292 /*@
4293   TSGetOptionsPrefix - Sets the prefix used for searching for all
4294   `TS` options in the database.
4295 
4296   Not Collective
4297 
4298   Input Parameter:
4299 . ts - The `TS` context
4300 
4301   Output Parameter:
4302 . prefix - A pointer to the prefix string used
4303 
4304   Level: intermediate
4305 
4306   Fortran Notes:
4307   The user should pass in a string 'prefix' of
4308   sufficient length to hold the prefix.
4309 
4310 .seealso: [](ch_ts), `TS`, `TSAppendOptionsPrefix()`, `TSSetFromOptions()`
4311 @*/
4312 PetscErrorCode TSGetOptionsPrefix(TS ts, const char *prefix[])
4313 {
4314   PetscFunctionBegin;
4315   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4316   PetscAssertPointer(prefix, 2);
4317   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ts, prefix));
4318   PetscFunctionReturn(PETSC_SUCCESS);
4319 }
4320 
4321 /*@C
4322   TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
4323 
4324   Not Collective, but parallel objects are returned if ts is parallel
4325 
4326   Input Parameter:
4327 . ts - The `TS` context obtained from `TSCreate()`
4328 
4329   Output Parameters:
4330 + Amat - The (approximate) Jacobian J of G, where U_t = G(U,t)  (or `NULL`)
4331 . Pmat - The matrix from which the preconditioner is constructed, usually the same as `Amat`  (or `NULL`)
4332 . func - Function to compute the Jacobian of the RHS  (or `NULL`)
4333 - ctx  - User-defined context for Jacobian evaluation routine  (or `NULL`)
4334 
4335   Level: intermediate
4336 
4337   Note:
4338   You can pass in `NULL` for any return argument you do not need.
4339 
4340 .seealso: [](ch_ts), `TS`, `TSGetTimeStep()`, `TSGetMatrices()`, `TSGetTime()`, `TSGetStepNumber()`
4341 
4342 @*/
4343 PetscErrorCode TSGetRHSJacobian(TS ts, Mat *Amat, Mat *Pmat, TSRHSJacobianFn **func, void **ctx)
4344 {
4345   DM dm;
4346 
4347   PetscFunctionBegin;
4348   if (Amat || Pmat) {
4349     SNES snes;
4350     PetscCall(TSGetSNES(ts, &snes));
4351     PetscCall(SNESSetUpMatrices(snes));
4352     PetscCall(SNESGetJacobian(snes, Amat, Pmat, NULL, NULL));
4353   }
4354   PetscCall(TSGetDM(ts, &dm));
4355   PetscCall(DMTSGetRHSJacobian(dm, func, ctx));
4356   PetscFunctionReturn(PETSC_SUCCESS);
4357 }
4358 
4359 /*@C
4360   TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
4361 
4362   Not Collective, but parallel objects are returned if ts is parallel
4363 
4364   Input Parameter:
4365 . ts - The `TS` context obtained from `TSCreate()`
4366 
4367   Output Parameters:
4368 + Amat - The (approximate) Jacobian of F(t,U,U_t)
4369 . Pmat - The matrix from which the preconditioner is constructed, often the same as `Amat`
4370 . f    - The function to compute the matrices
4371 - ctx  - User-defined context for Jacobian evaluation routine
4372 
4373   Level: advanced
4374 
4375   Note:
4376   You can pass in `NULL` for any return argument you do not need.
4377 
4378 .seealso: [](ch_ts), `TS`, `TSGetTimeStep()`, `TSGetRHSJacobian()`, `TSGetMatrices()`, `TSGetTime()`, `TSGetStepNumber()`
4379 @*/
4380 PetscErrorCode TSGetIJacobian(TS ts, Mat *Amat, Mat *Pmat, TSIJacobianFn **f, void **ctx)
4381 {
4382   DM dm;
4383 
4384   PetscFunctionBegin;
4385   if (Amat || Pmat) {
4386     SNES snes;
4387     PetscCall(TSGetSNES(ts, &snes));
4388     PetscCall(SNESSetUpMatrices(snes));
4389     PetscCall(SNESGetJacobian(snes, Amat, Pmat, NULL, NULL));
4390   }
4391   PetscCall(TSGetDM(ts, &dm));
4392   PetscCall(DMTSGetIJacobian(dm, f, ctx));
4393   PetscFunctionReturn(PETSC_SUCCESS);
4394 }
4395 
4396 #include <petsc/private/dmimpl.h>
4397 /*@
4398   TSSetDM - Sets the `DM` that may be used by some nonlinear solvers or preconditioners under the `TS`
4399 
4400   Logically Collective
4401 
4402   Input Parameters:
4403 + ts - the `TS` integrator object
4404 - dm - the dm, cannot be `NULL`
4405 
4406   Level: intermediate
4407 
4408   Notes:
4409   A `DM` can only be used for solving one problem at a time because information about the problem is stored on the `DM`,
4410   even when not using interfaces like `DMTSSetIFunction()`.  Use `DMClone()` to get a distinct `DM` when solving
4411   different problems using the same function space.
4412 
4413 .seealso: [](ch_ts), `TS`, `DM`, `TSGetDM()`, `SNESSetDM()`, `SNESGetDM()`
4414 @*/
4415 PetscErrorCode TSSetDM(TS ts, DM dm)
4416 {
4417   SNES snes;
4418   DMTS tsdm;
4419 
4420   PetscFunctionBegin;
4421   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4422   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
4423   PetscCall(PetscObjectReference((PetscObject)dm));
4424   if (ts->dm) { /* Move the DMTS context over to the new DM unless the new DM already has one */
4425     if (ts->dm->dmts && !dm->dmts) {
4426       PetscCall(DMCopyDMTS(ts->dm, dm));
4427       PetscCall(DMGetDMTS(ts->dm, &tsdm));
4428       /* Grant write privileges to the replacement DM */
4429       if (tsdm->originaldm == ts->dm) tsdm->originaldm = dm;
4430     }
4431     PetscCall(DMDestroy(&ts->dm));
4432   }
4433   ts->dm = dm;
4434 
4435   PetscCall(TSGetSNES(ts, &snes));
4436   PetscCall(SNESSetDM(snes, dm));
4437   PetscFunctionReturn(PETSC_SUCCESS);
4438 }
4439 
4440 /*@
4441   TSGetDM - Gets the `DM` that may be used by some preconditioners
4442 
4443   Not Collective
4444 
4445   Input Parameter:
4446 . ts - the `TS`
4447 
4448   Output Parameter:
4449 . dm - the `DM`
4450 
4451   Level: intermediate
4452 
4453 .seealso: [](ch_ts), `TS`, `DM`, `TSSetDM()`, `SNESSetDM()`, `SNESGetDM()`
4454 @*/
4455 PetscErrorCode TSGetDM(TS ts, DM *dm)
4456 {
4457   PetscFunctionBegin;
4458   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4459   if (!ts->dm) {
4460     PetscCall(DMShellCreate(PetscObjectComm((PetscObject)ts), &ts->dm));
4461     if (ts->snes) PetscCall(SNESSetDM(ts->snes, ts->dm));
4462   }
4463   *dm = ts->dm;
4464   PetscFunctionReturn(PETSC_SUCCESS);
4465 }
4466 
4467 /*@
4468   SNESTSFormFunction - Function to evaluate nonlinear residual
4469 
4470   Logically Collective
4471 
4472   Input Parameters:
4473 + snes - nonlinear solver
4474 . U    - the current state at which to evaluate the residual
4475 - ctx  - user context, must be a TS
4476 
4477   Output Parameter:
4478 . F - the nonlinear residual
4479 
4480   Level: advanced
4481 
4482   Note:
4483   This function is not normally called by users and is automatically registered with the `SNES` used by `TS`.
4484   It is most frequently passed to `MatFDColoringSetFunction()`.
4485 
4486 .seealso: [](ch_ts), `SNESSetFunction()`, `MatFDColoringSetFunction()`
4487 @*/
4488 PetscErrorCode SNESTSFormFunction(SNES snes, Vec U, Vec F, void *ctx)
4489 {
4490   TS ts = (TS)ctx;
4491 
4492   PetscFunctionBegin;
4493   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
4494   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
4495   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
4496   PetscValidHeaderSpecific(ts, TS_CLASSID, 4);
4497   PetscCheck(ts->ops->snesfunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "No method snesfunction for TS of type %s", ((PetscObject)ts)->type_name);
4498   PetscCall((*ts->ops->snesfunction)(snes, U, F, ts));
4499   PetscFunctionReturn(PETSC_SUCCESS);
4500 }
4501 
4502 /*@
4503   SNESTSFormJacobian - Function to evaluate the Jacobian
4504 
4505   Collective
4506 
4507   Input Parameters:
4508 + snes - nonlinear solver
4509 . U    - the current state at which to evaluate the residual
4510 - ctx  - user context, must be a `TS`
4511 
4512   Output Parameters:
4513 + A - the Jacobian
4514 - B - the preconditioning matrix (may be the same as A)
4515 
4516   Level: developer
4517 
4518   Note:
4519   This function is not normally called by users and is automatically registered with the `SNES` used by `TS`.
4520 
4521 .seealso: [](ch_ts), `SNESSetJacobian()`
4522 @*/
4523 PetscErrorCode SNESTSFormJacobian(SNES snes, Vec U, Mat A, Mat B, void *ctx)
4524 {
4525   TS ts = (TS)ctx;
4526 
4527   PetscFunctionBegin;
4528   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
4529   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
4530   PetscValidHeaderSpecific(A, MAT_CLASSID, 3);
4531   PetscValidHeaderSpecific(B, MAT_CLASSID, 4);
4532   PetscValidHeaderSpecific(ts, TS_CLASSID, 5);
4533   PetscCheck(ts->ops->snesjacobian, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "No method snesjacobian for TS of type %s", ((PetscObject)ts)->type_name);
4534   PetscCall((*ts->ops->snesjacobian)(snes, U, A, B, ts));
4535   PetscFunctionReturn(PETSC_SUCCESS);
4536 }
4537 
4538 /*@C
4539   TSComputeRHSFunctionLinear - Evaluate the right-hand side via the user-provided Jacobian, for linear problems Udot = A U only
4540 
4541   Collective
4542 
4543   Input Parameters:
4544 + ts  - time stepping context
4545 . t   - time at which to evaluate
4546 . U   - state at which to evaluate
4547 - ctx - context
4548 
4549   Output Parameter:
4550 . F - right-hand side
4551 
4552   Level: intermediate
4553 
4554   Note:
4555   This function is intended to be passed to `TSSetRHSFunction()` to evaluate the right-hand side for linear problems.
4556   The matrix (and optionally the evaluation context) should be passed to `TSSetRHSJacobian()`.
4557 
4558 .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `TSSetRHSJacobian()`, `TSComputeRHSJacobianConstant()`
4559 @*/
4560 PetscErrorCode TSComputeRHSFunctionLinear(TS ts, PetscReal t, Vec U, Vec F, void *ctx)
4561 {
4562   Mat Arhs, Brhs;
4563 
4564   PetscFunctionBegin;
4565   PetscCall(TSGetRHSMats_Private(ts, &Arhs, &Brhs));
4566   /* undo the damage caused by shifting */
4567   PetscCall(TSRecoverRHSJacobian(ts, Arhs, Brhs));
4568   PetscCall(TSComputeRHSJacobian(ts, t, U, Arhs, Brhs));
4569   PetscCall(MatMult(Arhs, U, F));
4570   PetscFunctionReturn(PETSC_SUCCESS);
4571 }
4572 
4573 /*@C
4574   TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
4575 
4576   Collective
4577 
4578   Input Parameters:
4579 + ts  - time stepping context
4580 . t   - time at which to evaluate
4581 . U   - state at which to evaluate
4582 - ctx - context
4583 
4584   Output Parameters:
4585 + A - pointer to operator
4586 - B - pointer to preconditioning matrix
4587 
4588   Level: intermediate
4589 
4590   Note:
4591   This function is intended to be passed to `TSSetRHSJacobian()` to evaluate the Jacobian for linear time-independent problems.
4592 
4593 .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `TSSetRHSJacobian()`, `TSComputeRHSFunctionLinear()`
4594 @*/
4595 PetscErrorCode TSComputeRHSJacobianConstant(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx)
4596 {
4597   PetscFunctionBegin;
4598   PetscFunctionReturn(PETSC_SUCCESS);
4599 }
4600 
4601 /*@C
4602   TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
4603 
4604   Collective
4605 
4606   Input Parameters:
4607 + ts   - time stepping context
4608 . t    - time at which to evaluate
4609 . U    - state at which to evaluate
4610 . Udot - time derivative of state vector
4611 - ctx  - context
4612 
4613   Output Parameter:
4614 . F - left hand side
4615 
4616   Level: intermediate
4617 
4618   Notes:
4619   The assumption here is that the left hand side is of the form A*Udot (and not A*Udot + B*U). For other cases, the
4620   user is required to write their own `TSComputeIFunction()`.
4621   This function is intended to be passed to `TSSetIFunction()` to evaluate the left hand side for linear problems.
4622   The matrix (and optionally the evaluation context) should be passed to `TSSetIJacobian()`.
4623 
4624   Note that using this function is NOT equivalent to using `TSComputeRHSFunctionLinear()` since that solves Udot = A U
4625 
4626 .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `TSSetIJacobian()`, `TSComputeIJacobianConstant()`, `TSComputeRHSFunctionLinear()`
4627 @*/
4628 PetscErrorCode TSComputeIFunctionLinear(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
4629 {
4630   Mat A, B;
4631 
4632   PetscFunctionBegin;
4633   PetscCall(TSGetIJacobian(ts, &A, &B, NULL, NULL));
4634   PetscCall(TSComputeIJacobian(ts, t, U, Udot, 1.0, A, B, PETSC_TRUE));
4635   PetscCall(MatMult(A, Udot, F));
4636   PetscFunctionReturn(PETSC_SUCCESS);
4637 }
4638 
4639 /*@C
4640   TSComputeIJacobianConstant - Reuses the matrix previously computed with the provided `TSIJacobianFn` for a semi-implicit DAE or ODE
4641 
4642   Collective
4643 
4644   Input Parameters:
4645 + ts    - time stepping context
4646 . t     - time at which to evaluate
4647 . U     - state at which to evaluate
4648 . Udot  - time derivative of state vector
4649 . shift - shift to apply
4650 - ctx   - context
4651 
4652   Output Parameters:
4653 + A - pointer to operator
4654 - B - pointer to matrix from which the preconditioner is built (often `A`)
4655 
4656   Level: advanced
4657 
4658   Notes:
4659   This function is intended to be passed to `TSSetIJacobian()` to evaluate the Jacobian for linear time-independent problems.
4660 
4661   It is only appropriate for problems of the form
4662 
4663   $$
4664   M \dot{U} = F(U,t)
4665   $$
4666 
4667   where M is constant and F is non-stiff.  The user must pass M to `TSSetIJacobian()`.  The current implementation only
4668   works with IMEX time integration methods such as `TSROSW` and `TSARKIMEX`, since there is no support for de-constructing
4669   an implicit operator of the form
4670 
4671   $$
4672   shift*M + J
4673   $$
4674 
4675   where J is the Jacobian of -F(U).  Support may be added in a future version of PETSc, but for now, the user must store
4676   a copy of M or reassemble it when requested.
4677 
4678 .seealso: [](ch_ts), `TS`, `TSROSW`, `TSARKIMEX`, `TSSetIFunction()`, `TSSetIJacobian()`, `TSComputeIFunctionLinear()`
4679 @*/
4680 PetscErrorCode TSComputeIJacobianConstant(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat A, Mat B, void *ctx)
4681 {
4682   PetscFunctionBegin;
4683   PetscCall(MatScale(A, shift / ts->ijacobian.shift));
4684   ts->ijacobian.shift = shift;
4685   PetscFunctionReturn(PETSC_SUCCESS);
4686 }
4687 
4688 /*@
4689   TSGetEquationType - Gets the type of the equation that `TS` is solving.
4690 
4691   Not Collective
4692 
4693   Input Parameter:
4694 . ts - the `TS` context
4695 
4696   Output Parameter:
4697 . equation_type - see `TSEquationType`
4698 
4699   Level: beginner
4700 
4701 .seealso: [](ch_ts), `TS`, `TSSetEquationType()`, `TSEquationType`
4702 @*/
4703 PetscErrorCode TSGetEquationType(TS ts, TSEquationType *equation_type)
4704 {
4705   PetscFunctionBegin;
4706   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4707   PetscAssertPointer(equation_type, 2);
4708   *equation_type = ts->equation_type;
4709   PetscFunctionReturn(PETSC_SUCCESS);
4710 }
4711 
4712 /*@
4713   TSSetEquationType - Sets the type of the equation that `TS` is solving.
4714 
4715   Not Collective
4716 
4717   Input Parameters:
4718 + ts            - the `TS` context
4719 - equation_type - see `TSEquationType`
4720 
4721   Level: advanced
4722 
4723 .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSEquationType`
4724 @*/
4725 PetscErrorCode TSSetEquationType(TS ts, TSEquationType equation_type)
4726 {
4727   PetscFunctionBegin;
4728   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4729   ts->equation_type = equation_type;
4730   PetscFunctionReturn(PETSC_SUCCESS);
4731 }
4732 
4733 /*@
4734   TSGetConvergedReason - Gets the reason the `TS` iteration was stopped.
4735 
4736   Not Collective
4737 
4738   Input Parameter:
4739 . ts - the `TS` context
4740 
4741   Output Parameter:
4742 . reason - negative value indicates diverged, positive value converged, see `TSConvergedReason` or the
4743             manual pages for the individual convergence tests for complete lists
4744 
4745   Level: beginner
4746 
4747   Note:
4748   Can only be called after the call to `TSSolve()` is complete.
4749 
4750 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSConvergedReason`
4751 @*/
4752 PetscErrorCode TSGetConvergedReason(TS ts, TSConvergedReason *reason)
4753 {
4754   PetscFunctionBegin;
4755   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4756   PetscAssertPointer(reason, 2);
4757   *reason = ts->reason;
4758   PetscFunctionReturn(PETSC_SUCCESS);
4759 }
4760 
4761 /*@
4762   TSSetConvergedReason - Sets the reason for handling the convergence of `TSSolve()`.
4763 
4764   Logically Collective; reason must contain common value
4765 
4766   Input Parameters:
4767 + ts     - the `TS` context
4768 - reason - negative value indicates diverged, positive value converged, see `TSConvergedReason` or the
4769             manual pages for the individual convergence tests for complete lists
4770 
4771   Level: advanced
4772 
4773   Note:
4774   Can only be called while `TSSolve()` is active.
4775 
4776 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSConvergedReason`
4777 @*/
4778 PetscErrorCode TSSetConvergedReason(TS ts, TSConvergedReason reason)
4779 {
4780   PetscFunctionBegin;
4781   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4782   ts->reason = reason;
4783   PetscFunctionReturn(PETSC_SUCCESS);
4784 }
4785 
4786 /*@
4787   TSGetSolveTime - Gets the time after a call to `TSSolve()`
4788 
4789   Not Collective
4790 
4791   Input Parameter:
4792 . ts - the `TS` context
4793 
4794   Output Parameter:
4795 . ftime - the final time. This time corresponds to the final time set with `TSSetMaxTime()`
4796 
4797   Level: beginner
4798 
4799   Note:
4800   Can only be called after the call to `TSSolve()` is complete.
4801 
4802 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSConvergedReason`
4803 @*/
4804 PetscErrorCode TSGetSolveTime(TS ts, PetscReal *ftime)
4805 {
4806   PetscFunctionBegin;
4807   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4808   PetscAssertPointer(ftime, 2);
4809   *ftime = ts->solvetime;
4810   PetscFunctionReturn(PETSC_SUCCESS);
4811 }
4812 
4813 /*@
4814   TSGetSNESIterations - Gets the total number of nonlinear iterations
4815   used by the time integrator.
4816 
4817   Not Collective
4818 
4819   Input Parameter:
4820 . ts - `TS` context
4821 
4822   Output Parameter:
4823 . nits - number of nonlinear iterations
4824 
4825   Level: intermediate
4826 
4827   Note:
4828   This counter is reset to zero for each successive call to `TSSolve()`.
4829 
4830 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetKSPIterations()`
4831 @*/
4832 PetscErrorCode TSGetSNESIterations(TS ts, PetscInt *nits)
4833 {
4834   PetscFunctionBegin;
4835   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4836   PetscAssertPointer(nits, 2);
4837   *nits = ts->snes_its;
4838   PetscFunctionReturn(PETSC_SUCCESS);
4839 }
4840 
4841 /*@
4842   TSGetKSPIterations - Gets the total number of linear iterations
4843   used by the time integrator.
4844 
4845   Not Collective
4846 
4847   Input Parameter:
4848 . ts - `TS` context
4849 
4850   Output Parameter:
4851 . lits - number of linear iterations
4852 
4853   Level: intermediate
4854 
4855   Note:
4856   This counter is reset to zero for each successive call to `TSSolve()`.
4857 
4858 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetSNESIterations()`, `SNESGetKSPIterations()`
4859 @*/
4860 PetscErrorCode TSGetKSPIterations(TS ts, PetscInt *lits)
4861 {
4862   PetscFunctionBegin;
4863   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4864   PetscAssertPointer(lits, 2);
4865   *lits = ts->ksp_its;
4866   PetscFunctionReturn(PETSC_SUCCESS);
4867 }
4868 
4869 /*@
4870   TSGetStepRejections - Gets the total number of rejected steps.
4871 
4872   Not Collective
4873 
4874   Input Parameter:
4875 . ts - `TS` context
4876 
4877   Output Parameter:
4878 . rejects - number of steps rejected
4879 
4880   Level: intermediate
4881 
4882   Note:
4883   This counter is reset to zero for each successive call to `TSSolve()`.
4884 
4885 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetSNESIterations()`, `TSGetKSPIterations()`, `TSSetMaxStepRejections()`, `TSGetSNESFailures()`, `TSSetMaxSNESFailures()`, `TSSetErrorIfStepFails()`
4886 @*/
4887 PetscErrorCode TSGetStepRejections(TS ts, PetscInt *rejects)
4888 {
4889   PetscFunctionBegin;
4890   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4891   PetscAssertPointer(rejects, 2);
4892   *rejects = ts->reject;
4893   PetscFunctionReturn(PETSC_SUCCESS);
4894 }
4895 
4896 /*@
4897   TSGetSNESFailures - Gets the total number of failed `SNES` solves in a `TS`
4898 
4899   Not Collective
4900 
4901   Input Parameter:
4902 . ts - `TS` context
4903 
4904   Output Parameter:
4905 . fails - number of failed nonlinear solves
4906 
4907   Level: intermediate
4908 
4909   Note:
4910   This counter is reset to zero for each successive call to `TSSolve()`.
4911 
4912 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetSNESIterations()`, `TSGetKSPIterations()`, `TSSetMaxStepRejections()`, `TSGetStepRejections()`, `TSSetMaxSNESFailures()`
4913 @*/
4914 PetscErrorCode TSGetSNESFailures(TS ts, PetscInt *fails)
4915 {
4916   PetscFunctionBegin;
4917   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4918   PetscAssertPointer(fails, 2);
4919   *fails = ts->num_snes_failures;
4920   PetscFunctionReturn(PETSC_SUCCESS);
4921 }
4922 
4923 /*@
4924   TSSetMaxStepRejections - Sets the maximum number of step rejections before a time step fails
4925 
4926   Not Collective
4927 
4928   Input Parameters:
4929 + ts      - `TS` context
4930 - rejects - maximum number of rejected steps, pass `PETSC_UNLIMITED` for unlimited
4931 
4932   Options Database Key:
4933 . -ts_max_reject - Maximum number of step rejections before a step fails
4934 
4935   Level: intermediate
4936 
4937   Developer Note:
4938   The options database name is incorrect.
4939 
4940 .seealso: [](ch_ts), `TS`, `SNES`, `TSGetSNESIterations()`, `TSGetKSPIterations()`, `TSSetMaxSNESFailures()`, `TSGetStepRejections()`, `TSGetSNESFailures()`, `TSSetErrorIfStepFails()`, `TSGetConvergedReason()`
4941 @*/
4942 PetscErrorCode TSSetMaxStepRejections(TS ts, PetscInt rejects)
4943 {
4944   PetscFunctionBegin;
4945   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4946   if (rejects == PETSC_UNLIMITED || rejects == -1) {
4947     ts->max_reject = PETSC_UNLIMITED;
4948   } else {
4949     PetscCheck(rejects >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have a negative maximum number of rejections");
4950     ts->max_reject = rejects;
4951   }
4952   PetscFunctionReturn(PETSC_SUCCESS);
4953 }
4954 
4955 /*@
4956   TSSetMaxSNESFailures - Sets the maximum number of failed `SNES` solves
4957 
4958   Not Collective
4959 
4960   Input Parameters:
4961 + ts    - `TS` context
4962 - fails - maximum number of failed nonlinear solves, pass `PETSC_UNLIMITED` to allow any number of failures.
4963 
4964   Options Database Key:
4965 . -ts_max_snes_failures - Maximum number of nonlinear solve failures
4966 
4967   Level: intermediate
4968 
4969 .seealso: [](ch_ts), `TS`, `SNES`, `TSGetSNESIterations()`, `TSGetKSPIterations()`, `TSSetMaxStepRejections()`, `TSGetStepRejections()`, `TSGetSNESFailures()`, `SNESGetConvergedReason()`, `TSGetConvergedReason()`
4970 @*/
4971 PetscErrorCode TSSetMaxSNESFailures(TS ts, PetscInt fails)
4972 {
4973   PetscFunctionBegin;
4974   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4975   if (fails == PETSC_UNLIMITED || fails == -1) {
4976     ts->max_snes_failures = PETSC_UNLIMITED;
4977   } else {
4978     PetscCheck(fails >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have a negative maximum number of failures");
4979     ts->max_snes_failures = fails;
4980   }
4981   PetscFunctionReturn(PETSC_SUCCESS);
4982 }
4983 
4984 /*@
4985   TSSetErrorIfStepFails - Immediately error if no step succeeds during `TSSolve()`
4986 
4987   Not Collective
4988 
4989   Input Parameters:
4990 + ts  - `TS` context
4991 - err - `PETSC_TRUE` to error if no step succeeds, `PETSC_FALSE` to return without failure
4992 
4993   Options Database Key:
4994 . -ts_error_if_step_fails - Error if no step succeeds
4995 
4996   Level: intermediate
4997 
4998 .seealso: [](ch_ts), `TS`, `TSGetSNESIterations()`, `TSGetKSPIterations()`, `TSSetMaxStepRejections()`, `TSGetStepRejections()`, `TSGetSNESFailures()`, `TSGetConvergedReason()`
4999 @*/
5000 PetscErrorCode TSSetErrorIfStepFails(TS ts, PetscBool err)
5001 {
5002   PetscFunctionBegin;
5003   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5004   ts->errorifstepfailed = err;
5005   PetscFunctionReturn(PETSC_SUCCESS);
5006 }
5007 
5008 /*@
5009   TSGetAdapt - Get the adaptive controller context for the current method
5010 
5011   Collective if controller has not yet been created
5012 
5013   Input Parameter:
5014 . ts - time stepping context
5015 
5016   Output Parameter:
5017 . adapt - adaptive controller
5018 
5019   Level: intermediate
5020 
5021 .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSAdaptSetType()`, `TSAdaptChoose()`
5022 @*/
5023 PetscErrorCode TSGetAdapt(TS ts, TSAdapt *adapt)
5024 {
5025   PetscFunctionBegin;
5026   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5027   PetscAssertPointer(adapt, 2);
5028   if (!ts->adapt) {
5029     PetscCall(TSAdaptCreate(PetscObjectComm((PetscObject)ts), &ts->adapt));
5030     PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->adapt, (PetscObject)ts, 1));
5031   }
5032   *adapt = ts->adapt;
5033   PetscFunctionReturn(PETSC_SUCCESS);
5034 }
5035 
5036 /*@
5037   TSSetTolerances - Set tolerances for local truncation error when using an adaptive controller
5038 
5039   Logically Collective
5040 
5041   Input Parameters:
5042 + ts    - time integration context
5043 . atol  - scalar absolute tolerances
5044 . vatol - vector of absolute tolerances or `NULL`, used in preference to `atol` if present
5045 . rtol  - scalar relative tolerances
5046 - vrtol - vector of relative tolerances or `NULL`, used in preference to `rtol` if present
5047 
5048   Options Database Keys:
5049 + -ts_rtol <rtol> - relative tolerance for local truncation error
5050 - -ts_atol <atol> - Absolute tolerance for local truncation error
5051 
5052   Level: beginner
5053 
5054   Notes:
5055   `PETSC_CURRENT` or `PETSC_DETERMINE` may be used for `atol` or `rtol` to indicate the current value
5056   or the default value from when the object's type was set.
5057 
5058   With PETSc's implicit schemes for DAE problems, the calculation of the local truncation error
5059   (LTE) includes both the differential and the algebraic variables. If one wants the LTE to be
5060   computed only for the differential or the algebraic part then this can be done using the vector of
5061   tolerances vatol. For example, by setting the tolerance vector with the desired tolerance for the
5062   differential part and infinity for the algebraic part, the LTE calculation will include only the
5063   differential variables.
5064 
5065   Fortran Note:
5066   Use `PETSC_CURRENT_INTEGER` or `PETSC_DETERMINE_INTEGER`.
5067 
5068 .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSErrorWeightedNorm()`, `TSGetTolerances()`
5069 @*/
5070 PetscErrorCode TSSetTolerances(TS ts, PetscReal atol, Vec vatol, PetscReal rtol, Vec vrtol)
5071 {
5072   PetscFunctionBegin;
5073   if (atol == (PetscReal)PETSC_DETERMINE) {
5074     ts->atol = ts->default_atol;
5075   } else if (atol != (PetscReal)PETSC_CURRENT) {
5076     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)atol);
5077     ts->atol = atol;
5078   }
5079 
5080   if (vatol) {
5081     PetscCall(PetscObjectReference((PetscObject)vatol));
5082     PetscCall(VecDestroy(&ts->vatol));
5083     ts->vatol = vatol;
5084   }
5085 
5086   if (rtol == (PetscReal)PETSC_DETERMINE) {
5087     ts->rtol = ts->default_rtol;
5088   } else if (rtol != (PetscReal)PETSC_CURRENT) {
5089     PetscCheck(rtol >= 0.0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Relative tolerance %g must be non-negative", (double)rtol);
5090     ts->rtol = rtol;
5091   }
5092 
5093   if (vrtol) {
5094     PetscCall(PetscObjectReference((PetscObject)vrtol));
5095     PetscCall(VecDestroy(&ts->vrtol));
5096     ts->vrtol = vrtol;
5097   }
5098   PetscFunctionReturn(PETSC_SUCCESS);
5099 }
5100 
5101 /*@
5102   TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
5103 
5104   Logically Collective
5105 
5106   Input Parameter:
5107 . ts - time integration context
5108 
5109   Output Parameters:
5110 + atol  - scalar absolute tolerances, `NULL` to ignore
5111 . vatol - vector of absolute tolerances, `NULL` to ignore
5112 . rtol  - scalar relative tolerances, `NULL` to ignore
5113 - vrtol - vector of relative tolerances, `NULL` to ignore
5114 
5115   Level: beginner
5116 
5117 .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSErrorWeightedNorm()`, `TSSetTolerances()`
5118 @*/
5119 PetscErrorCode TSGetTolerances(TS ts, PetscReal *atol, Vec *vatol, PetscReal *rtol, Vec *vrtol)
5120 {
5121   PetscFunctionBegin;
5122   if (atol) *atol = ts->atol;
5123   if (vatol) *vatol = ts->vatol;
5124   if (rtol) *rtol = ts->rtol;
5125   if (vrtol) *vrtol = ts->vrtol;
5126   PetscFunctionReturn(PETSC_SUCCESS);
5127 }
5128 
5129 /*@
5130   TSErrorWeightedNorm - compute a weighted norm of the difference between two state vectors based on supplied absolute and relative tolerances
5131 
5132   Collective
5133 
5134   Input Parameters:
5135 + ts        - time stepping context
5136 . U         - state vector, usually ts->vec_sol
5137 . Y         - state vector to be compared to U
5138 - wnormtype - norm type, either `NORM_2` or `NORM_INFINITY`
5139 
5140   Output Parameters:
5141 + norm  - weighted norm, a value of 1.0 achieves a balance between absolute and relative tolerances
5142 . norma - weighted norm, a value of 1.0 means that the error meets the absolute tolerance set by the user
5143 - normr - weighted norm, a value of 1.0 means that the error meets the relative tolerance set by the user
5144 
5145   Options Database Key:
5146 . -ts_adapt_wnormtype <wnormtype> - 2, INFINITY
5147 
5148   Level: developer
5149 
5150 .seealso: [](ch_ts), `TS`, `VecErrorWeightedNorms()`, `TSErrorWeightedENorm()`
5151 @*/
5152 PetscErrorCode TSErrorWeightedNorm(TS ts, Vec U, Vec Y, NormType wnormtype, PetscReal *norm, PetscReal *norma, PetscReal *normr)
5153 {
5154   PetscInt norma_loc, norm_loc, normr_loc;
5155 
5156   PetscFunctionBegin;
5157   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5158   PetscCall(VecErrorWeightedNorms(U, Y, NULL, wnormtype, ts->atol, ts->vatol, ts->rtol, ts->vrtol, ts->adapt->ignore_max, norm, &norm_loc, norma, &norma_loc, normr, &normr_loc));
5159   if (wnormtype == NORM_2) {
5160     if (norm_loc) *norm = PetscSqrtReal(PetscSqr(*norm) / norm_loc);
5161     if (norma_loc) *norma = PetscSqrtReal(PetscSqr(*norma) / norma_loc);
5162     if (normr_loc) *normr = PetscSqrtReal(PetscSqr(*normr) / normr_loc);
5163   }
5164   PetscCheck(!PetscIsInfOrNanScalar(*norm), PetscObjectComm((PetscObject)ts), PETSC_ERR_FP, "Infinite or not-a-number generated in norm");
5165   PetscCheck(!PetscIsInfOrNanScalar(*norma), PetscObjectComm((PetscObject)ts), PETSC_ERR_FP, "Infinite or not-a-number generated in norma");
5166   PetscCheck(!PetscIsInfOrNanScalar(*normr), PetscObjectComm((PetscObject)ts), PETSC_ERR_FP, "Infinite or not-a-number generated in normr");
5167   PetscFunctionReturn(PETSC_SUCCESS);
5168 }
5169 
5170 /*@
5171   TSErrorWeightedENorm - compute a weighted error norm based on supplied absolute and relative tolerances
5172 
5173   Collective
5174 
5175   Input Parameters:
5176 + ts        - time stepping context
5177 . E         - error vector
5178 . U         - state vector, usually ts->vec_sol
5179 . Y         - state vector, previous time step
5180 - wnormtype - norm type, either `NORM_2` or `NORM_INFINITY`
5181 
5182   Output Parameters:
5183 + norm  - weighted norm, a value of 1.0 achieves a balance between absolute and relative tolerances
5184 . norma - weighted norm, a value of 1.0 means that the error meets the absolute tolerance set by the user
5185 - normr - weighted norm, a value of 1.0 means that the error meets the relative tolerance set by the user
5186 
5187   Options Database Key:
5188 . -ts_adapt_wnormtype <wnormtype> - 2, INFINITY
5189 
5190   Level: developer
5191 
5192 .seealso: [](ch_ts), `TS`, `VecErrorWeightedNorms()`, `TSErrorWeightedNorm()`
5193 @*/
5194 PetscErrorCode TSErrorWeightedENorm(TS ts, Vec E, Vec U, Vec Y, NormType wnormtype, PetscReal *norm, PetscReal *norma, PetscReal *normr)
5195 {
5196   PetscInt norma_loc, norm_loc, normr_loc;
5197 
5198   PetscFunctionBegin;
5199   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5200   PetscCall(VecErrorWeightedNorms(U, Y, E, wnormtype, ts->atol, ts->vatol, ts->rtol, ts->vrtol, ts->adapt->ignore_max, norm, &norm_loc, norma, &norma_loc, normr, &normr_loc));
5201   if (wnormtype == NORM_2) {
5202     if (norm_loc) *norm = PetscSqrtReal(PetscSqr(*norm) / norm_loc);
5203     if (norma_loc) *norma = PetscSqrtReal(PetscSqr(*norma) / norma_loc);
5204     if (normr_loc) *normr = PetscSqrtReal(PetscSqr(*normr) / normr_loc);
5205   }
5206   PetscCheck(!PetscIsInfOrNanScalar(*norm), PetscObjectComm((PetscObject)ts), PETSC_ERR_FP, "Infinite or not-a-number generated in norm");
5207   PetscCheck(!PetscIsInfOrNanScalar(*norma), PetscObjectComm((PetscObject)ts), PETSC_ERR_FP, "Infinite or not-a-number generated in norma");
5208   PetscCheck(!PetscIsInfOrNanScalar(*normr), PetscObjectComm((PetscObject)ts), PETSC_ERR_FP, "Infinite or not-a-number generated in normr");
5209   PetscFunctionReturn(PETSC_SUCCESS);
5210 }
5211 
5212 /*@
5213   TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
5214 
5215   Logically Collective
5216 
5217   Input Parameters:
5218 + ts      - time stepping context
5219 - cfltime - maximum stable time step if using forward Euler (value can be different on each process)
5220 
5221   Note:
5222   After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
5223 
5224   Level: intermediate
5225 
5226 .seealso: [](ch_ts), `TSGetCFLTime()`, `TSADAPTCFL`
5227 @*/
5228 PetscErrorCode TSSetCFLTimeLocal(TS ts, PetscReal cfltime)
5229 {
5230   PetscFunctionBegin;
5231   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5232   ts->cfltime_local = cfltime;
5233   ts->cfltime       = -1.;
5234   PetscFunctionReturn(PETSC_SUCCESS);
5235 }
5236 
5237 /*@
5238   TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
5239 
5240   Collective
5241 
5242   Input Parameter:
5243 . ts - time stepping context
5244 
5245   Output Parameter:
5246 . cfltime - maximum stable time step for forward Euler
5247 
5248   Level: advanced
5249 
5250 .seealso: [](ch_ts), `TSSetCFLTimeLocal()`
5251 @*/
5252 PetscErrorCode TSGetCFLTime(TS ts, PetscReal *cfltime)
5253 {
5254   PetscFunctionBegin;
5255   if (ts->cfltime < 0) PetscCallMPI(MPIU_Allreduce(&ts->cfltime_local, &ts->cfltime, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts)));
5256   *cfltime = ts->cfltime;
5257   PetscFunctionReturn(PETSC_SUCCESS);
5258 }
5259 
5260 /*@
5261   TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
5262 
5263   Input Parameters:
5264 + ts - the `TS` context.
5265 . xl - lower bound.
5266 - xu - upper bound.
5267 
5268   Level: advanced
5269 
5270   Note:
5271   If this routine is not called then the lower and upper bounds are set to
5272   `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`.
5273 
5274 .seealso: [](ch_ts), `TS`
5275 @*/
5276 PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
5277 {
5278   SNES snes;
5279 
5280   PetscFunctionBegin;
5281   PetscCall(TSGetSNES(ts, &snes));
5282   PetscCall(SNESVISetVariableBounds(snes, xl, xu));
5283   PetscFunctionReturn(PETSC_SUCCESS);
5284 }
5285 
5286 /*@
5287   TSComputeLinearStability - computes the linear stability function at a point
5288 
5289   Collective
5290 
5291   Input Parameters:
5292 + ts - the `TS` context
5293 . xr - real part of input argument
5294 - xi - imaginary part of input argument
5295 
5296   Output Parameters:
5297 + yr - real part of function value
5298 - yi - imaginary part of function value
5299 
5300   Level: developer
5301 
5302 .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `TSComputeIFunction()`
5303 @*/
5304 PetscErrorCode TSComputeLinearStability(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
5305 {
5306   PetscFunctionBegin;
5307   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5308   PetscUseTypeMethod(ts, linearstability, xr, xi, yr, yi);
5309   PetscFunctionReturn(PETSC_SUCCESS);
5310 }
5311 
5312 /*@
5313   TSRestartStep - Flags the solver to restart the next step
5314 
5315   Collective
5316 
5317   Input Parameter:
5318 . ts - the `TS` context obtained from `TSCreate()`
5319 
5320   Level: advanced
5321 
5322   Notes:
5323   Multistep methods like `TSBDF` or Runge-Kutta methods with FSAL property require restarting the solver in the event of
5324   discontinuities. These discontinuities may be introduced as a consequence of explicitly modifications to the solution
5325   vector (which PETSc attempts to detect and handle) or problem coefficients (which PETSc is not able to detect). For
5326   the sake of correctness and maximum safety, users are expected to call `TSRestart()` whenever they introduce
5327   discontinuities in callback routines (e.g. prestep and poststep routines, or implicit/rhs function routines with
5328   discontinuous source terms).
5329 
5330 .seealso: [](ch_ts), `TS`, `TSBDF`, `TSSolve()`, `TSSetPreStep()`, `TSSetPostStep()`
5331 @*/
5332 PetscErrorCode TSRestartStep(TS ts)
5333 {
5334   PetscFunctionBegin;
5335   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5336   ts->steprestart = PETSC_TRUE;
5337   PetscFunctionReturn(PETSC_SUCCESS);
5338 }
5339 
5340 /*@
5341   TSRollBack - Rolls back one time step
5342 
5343   Collective
5344 
5345   Input Parameter:
5346 . ts - the `TS` context obtained from `TSCreate()`
5347 
5348   Level: advanced
5349 
5350 .seealso: [](ch_ts), `TS`, `TSGetStepRollBack()`, `TSCreate()`, `TSSetUp()`, `TSDestroy()`, `TSSolve()`, `TSSetPreStep()`, `TSSetPreStage()`, `TSInterpolate()`
5351 @*/
5352 PetscErrorCode TSRollBack(TS ts)
5353 {
5354   PetscFunctionBegin;
5355   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5356   PetscCheck(!ts->steprollback, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "TSRollBack already called");
5357   PetscTryTypeMethod(ts, rollback);
5358   PetscCall(VecCopy(ts->vec_sol0, ts->vec_sol));
5359   ts->time_step  = ts->ptime - ts->ptime_prev;
5360   ts->ptime      = ts->ptime_prev;
5361   ts->ptime_prev = ts->ptime_prev_rollback;
5362   ts->steps--;
5363   ts->steprollback = PETSC_TRUE;
5364   PetscFunctionReturn(PETSC_SUCCESS);
5365 }
5366 
5367 /*@
5368   TSGetStepRollBack - Get the internal flag indicating if you are rolling back a step
5369 
5370   Not collective
5371 
5372   Input Parameter:
5373 . ts - the `TS` context obtained from `TSCreate()`
5374 
5375   Output Parameter:
5376 . flg - the rollback flag
5377 
5378   Level: advanced
5379 
5380 .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSRollBack()`
5381 @*/
5382 PetscErrorCode TSGetStepRollBack(TS ts, PetscBool *flg)
5383 {
5384   PetscFunctionBegin;
5385   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5386   PetscAssertPointer(flg, 2);
5387   *flg = ts->steprollback;
5388   PetscFunctionReturn(PETSC_SUCCESS);
5389 }
5390 
5391 /*@
5392   TSGetStepResize - Get the internal flag indicating if the current step is after a resize.
5393 
5394   Not collective
5395 
5396   Input Parameter:
5397 . ts - the `TS` context obtained from `TSCreate()`
5398 
5399   Output Parameter:
5400 . flg - the resize flag
5401 
5402   Level: advanced
5403 
5404 .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetResize()`
5405 @*/
5406 PetscErrorCode TSGetStepResize(TS ts, PetscBool *flg)
5407 {
5408   PetscFunctionBegin;
5409   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5410   PetscAssertPointer(flg, 2);
5411   *flg = ts->stepresize;
5412   PetscFunctionReturn(PETSC_SUCCESS);
5413 }
5414 
5415 /*@
5416   TSGetStages - Get the number of stages and stage values
5417 
5418   Input Parameter:
5419 . ts - the `TS` context obtained from `TSCreate()`
5420 
5421   Output Parameters:
5422 + ns - the number of stages
5423 - Y  - the current stage vectors
5424 
5425   Level: advanced
5426 
5427   Note:
5428   Both `ns` and `Y` can be `NULL`.
5429 
5430 .seealso: [](ch_ts), `TS`, `TSCreate()`
5431 @*/
5432 PetscErrorCode TSGetStages(TS ts, PetscInt *ns, Vec **Y)
5433 {
5434   PetscFunctionBegin;
5435   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5436   if (ns) PetscAssertPointer(ns, 2);
5437   if (Y) PetscAssertPointer(Y, 3);
5438   if (!ts->ops->getstages) {
5439     if (ns) *ns = 0;
5440     if (Y) *Y = NULL;
5441   } else PetscUseTypeMethod(ts, getstages, ns, Y);
5442   PetscFunctionReturn(PETSC_SUCCESS);
5443 }
5444 
5445 /*@C
5446   TSComputeIJacobianDefaultColor - Computes the Jacobian using finite differences and coloring to exploit matrix sparsity.
5447 
5448   Collective
5449 
5450   Input Parameters:
5451 + ts    - the `TS` context
5452 . t     - current timestep
5453 . U     - state vector
5454 . Udot  - time derivative of state vector
5455 . shift - shift to apply, see note below
5456 - ctx   - an optional user context
5457 
5458   Output Parameters:
5459 + J - Jacobian matrix (not altered in this routine)
5460 - B - newly computed Jacobian matrix to use with preconditioner (generally the same as `J`)
5461 
5462   Level: intermediate
5463 
5464   Notes:
5465   If F(t,U,Udot)=0 is the DAE, the required Jacobian is
5466 
5467   dF/dU + shift*dF/dUdot
5468 
5469   Most users should not need to explicitly call this routine, as it
5470   is used internally within the nonlinear solvers.
5471 
5472   This will first try to get the coloring from the `DM`.  If the `DM` type has no coloring
5473   routine, then it will try to get the coloring from the matrix.  This requires that the
5474   matrix have nonzero entries precomputed.
5475 
5476 .seealso: [](ch_ts), `TS`, `TSSetIJacobian()`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`
5477 @*/
5478 PetscErrorCode TSComputeIJacobianDefaultColor(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat J, Mat B, void *ctx)
5479 {
5480   SNES          snes;
5481   MatFDColoring color;
5482   PetscBool     hascolor, matcolor = PETSC_FALSE;
5483 
5484   PetscFunctionBegin;
5485   PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_fd_color_use_mat", &matcolor, NULL));
5486   PetscCall(PetscObjectQuery((PetscObject)B, "TSMatFDColoring", (PetscObject *)&color));
5487   if (!color) {
5488     DM         dm;
5489     ISColoring iscoloring;
5490 
5491     PetscCall(TSGetDM(ts, &dm));
5492     PetscCall(DMHasColoring(dm, &hascolor));
5493     if (hascolor && !matcolor) {
5494       PetscCall(DMCreateColoring(dm, IS_COLORING_GLOBAL, &iscoloring));
5495       PetscCall(MatFDColoringCreate(B, iscoloring, &color));
5496       PetscCall(MatFDColoringSetFunction(color, (PetscErrorCode (*)(void))SNESTSFormFunction, (void *)ts));
5497       PetscCall(MatFDColoringSetFromOptions(color));
5498       PetscCall(MatFDColoringSetUp(B, iscoloring, color));
5499       PetscCall(ISColoringDestroy(&iscoloring));
5500     } else {
5501       MatColoring mc;
5502 
5503       PetscCall(MatColoringCreate(B, &mc));
5504       PetscCall(MatColoringSetDistance(mc, 2));
5505       PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
5506       PetscCall(MatColoringSetFromOptions(mc));
5507       PetscCall(MatColoringApply(mc, &iscoloring));
5508       PetscCall(MatColoringDestroy(&mc));
5509       PetscCall(MatFDColoringCreate(B, iscoloring, &color));
5510       PetscCall(MatFDColoringSetFunction(color, (PetscErrorCode (*)(void))SNESTSFormFunction, (void *)ts));
5511       PetscCall(MatFDColoringSetFromOptions(color));
5512       PetscCall(MatFDColoringSetUp(B, iscoloring, color));
5513       PetscCall(ISColoringDestroy(&iscoloring));
5514     }
5515     PetscCall(PetscObjectCompose((PetscObject)B, "TSMatFDColoring", (PetscObject)color));
5516     PetscCall(PetscObjectDereference((PetscObject)color));
5517   }
5518   PetscCall(TSGetSNES(ts, &snes));
5519   PetscCall(MatFDColoringApply(B, color, U, snes));
5520   if (J != B) {
5521     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
5522     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
5523   }
5524   PetscFunctionReturn(PETSC_SUCCESS);
5525 }
5526 
5527 /*@C
5528   TSSetFunctionDomainError - Set a function that tests if the current state vector is valid
5529 
5530   Input Parameters:
5531 + ts   - the `TS` context
5532 - func - function called within `TSFunctionDomainError()`
5533 
5534   Calling sequence of `func`:
5535 + ts     - the `TS` context
5536 . time   - the current time (of the stage)
5537 . state  - the state to check if it is valid
5538 - accept - (output parameter) `PETSC_FALSE` if the state is not acceptable, `PETSC_TRUE` if acceptable
5539 
5540   Level: intermediate
5541 
5542   Notes:
5543   If an implicit ODE solver is being used then, in addition to providing this routine, the
5544   user's code should call `SNESSetFunctionDomainError()` when domain errors occur during
5545   function evaluations where the functions are provided by `TSSetIFunction()` or `TSSetRHSFunction()`.
5546   Use `TSGetSNES()` to obtain the `SNES` object
5547 
5548   Developer Notes:
5549   The naming of this function is inconsistent with the `SNESSetFunctionDomainError()`
5550   since one takes a function pointer and the other does not.
5551 
5552 .seealso: [](ch_ts), `TSAdaptCheckStage()`, `TSFunctionDomainError()`, `SNESSetFunctionDomainError()`, `TSGetSNES()`
5553 @*/
5554 PetscErrorCode TSSetFunctionDomainError(TS ts, PetscErrorCode (*func)(TS ts, PetscReal time, Vec state, PetscBool *accept))
5555 {
5556   PetscFunctionBegin;
5557   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5558   ts->functiondomainerror = func;
5559   PetscFunctionReturn(PETSC_SUCCESS);
5560 }
5561 
5562 /*@
5563   TSFunctionDomainError - Checks if the current state is valid
5564 
5565   Input Parameters:
5566 + ts        - the `TS` context
5567 . stagetime - time of the simulation
5568 - Y         - state vector to check.
5569 
5570   Output Parameter:
5571 . accept - Set to `PETSC_FALSE` if the current state vector is valid.
5572 
5573   Level: developer
5574 
5575   Note:
5576   This function is called by the `TS` integration routines and calls the user provided function (set with `TSSetFunctionDomainError()`)
5577   to check if the current state is valid.
5578 
5579 .seealso: [](ch_ts), `TS`, `TSSetFunctionDomainError()`
5580 @*/
5581 PetscErrorCode TSFunctionDomainError(TS ts, PetscReal stagetime, Vec Y, PetscBool *accept)
5582 {
5583   PetscFunctionBegin;
5584   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5585   *accept = PETSC_TRUE;
5586   if (ts->functiondomainerror) PetscCall((*ts->functiondomainerror)(ts, stagetime, Y, accept));
5587   PetscFunctionReturn(PETSC_SUCCESS);
5588 }
5589 
5590 /*@
5591   TSClone - This function clones a time step `TS` object.
5592 
5593   Collective
5594 
5595   Input Parameter:
5596 . tsin - The input `TS`
5597 
5598   Output Parameter:
5599 . tsout - The output `TS` (cloned)
5600 
5601   Level: developer
5602 
5603   Notes:
5604   This function is used to create a clone of a `TS` object. It is used in `TSARKIMEX` for initializing the slope for first stage explicit methods.
5605   It will likely be replaced in the future with a mechanism of switching methods on the fly.
5606 
5607   When using `TSDestroy()` on a clone the user has to first reset the correct `TS` reference in the embedded `SNES` object: e.g., by running
5608 .vb
5609  SNES snes_dup = NULL;
5610  TSGetSNES(ts,&snes_dup);
5611  TSSetSNES(ts,snes_dup);
5612 .ve
5613 
5614 .seealso: [](ch_ts), `TS`, `SNES`, `TSCreate()`, `TSSetType()`, `TSSetUp()`, `TSDestroy()`, `TSSetProblemType()`
5615 @*/
5616 PetscErrorCode TSClone(TS tsin, TS *tsout)
5617 {
5618   TS     t;
5619   SNES   snes_start;
5620   DM     dm;
5621   TSType type;
5622 
5623   PetscFunctionBegin;
5624   PetscAssertPointer(tsin, 1);
5625   *tsout = NULL;
5626 
5627   PetscCall(PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", PetscObjectComm((PetscObject)tsin), TSDestroy, TSView));
5628 
5629   /* General TS description */
5630   t->numbermonitors    = 0;
5631   t->monitorFrequency  = 1;
5632   t->setupcalled       = 0;
5633   t->ksp_its           = 0;
5634   t->snes_its          = 0;
5635   t->nwork             = 0;
5636   t->rhsjacobian.time  = PETSC_MIN_REAL;
5637   t->rhsjacobian.scale = 1.;
5638   t->ijacobian.shift   = 1.;
5639 
5640   PetscCall(TSGetSNES(tsin, &snes_start));
5641   PetscCall(TSSetSNES(t, snes_start));
5642 
5643   PetscCall(TSGetDM(tsin, &dm));
5644   PetscCall(TSSetDM(t, dm));
5645 
5646   t->adapt = tsin->adapt;
5647   PetscCall(PetscObjectReference((PetscObject)t->adapt));
5648 
5649   t->trajectory = tsin->trajectory;
5650   PetscCall(PetscObjectReference((PetscObject)t->trajectory));
5651 
5652   t->event = tsin->event;
5653   if (t->event) t->event->refct++;
5654 
5655   t->problem_type      = tsin->problem_type;
5656   t->ptime             = tsin->ptime;
5657   t->ptime_prev        = tsin->ptime_prev;
5658   t->time_step         = tsin->time_step;
5659   t->max_time          = tsin->max_time;
5660   t->steps             = tsin->steps;
5661   t->max_steps         = tsin->max_steps;
5662   t->equation_type     = tsin->equation_type;
5663   t->atol              = tsin->atol;
5664   t->rtol              = tsin->rtol;
5665   t->max_snes_failures = tsin->max_snes_failures;
5666   t->max_reject        = tsin->max_reject;
5667   t->errorifstepfailed = tsin->errorifstepfailed;
5668 
5669   PetscCall(TSGetType(tsin, &type));
5670   PetscCall(TSSetType(t, type));
5671 
5672   t->vec_sol = NULL;
5673 
5674   t->cfltime          = tsin->cfltime;
5675   t->cfltime_local    = tsin->cfltime_local;
5676   t->exact_final_time = tsin->exact_final_time;
5677 
5678   t->ops[0] = tsin->ops[0];
5679 
5680   if (((PetscObject)tsin)->fortran_func_pointers) {
5681     PetscInt i;
5682     PetscCall(PetscMalloc((10) * sizeof(void (*)(void)), &((PetscObject)t)->fortran_func_pointers));
5683     for (i = 0; i < 10; i++) ((PetscObject)t)->fortran_func_pointers[i] = ((PetscObject)tsin)->fortran_func_pointers[i];
5684   }
5685   *tsout = t;
5686   PetscFunctionReturn(PETSC_SUCCESS);
5687 }
5688 
5689 static PetscErrorCode RHSWrapperFunction_TSRHSJacobianTest(void *ctx, Vec x, Vec y)
5690 {
5691   TS ts = (TS)ctx;
5692 
5693   PetscFunctionBegin;
5694   PetscCall(TSComputeRHSFunction(ts, 0, x, y));
5695   PetscFunctionReturn(PETSC_SUCCESS);
5696 }
5697 
5698 /*@
5699   TSRHSJacobianTest - Compares the multiply routine provided to the `MATSHELL` with differencing on the `TS` given RHS function.
5700 
5701   Logically Collective
5702 
5703   Input Parameter:
5704 . ts - the time stepping routine
5705 
5706   Output Parameter:
5707 . flg - `PETSC_TRUE` if the multiply is likely correct
5708 
5709   Options Database Key:
5710 . -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view - run the test at each timestep of the integrator
5711 
5712   Level: advanced
5713 
5714   Note:
5715   This only works for problems defined using `TSSetRHSFunction()` and Jacobian NOT `TSSetIFunction()` and Jacobian
5716 
5717 .seealso: [](ch_ts), `TS`, `Mat`, `MATSHELL`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellTestMultTranspose()`, `TSRHSJacobianTestTranspose()`
5718 @*/
5719 PetscErrorCode TSRHSJacobianTest(TS ts, PetscBool *flg)
5720 {
5721   Mat              J, B;
5722   TSRHSJacobianFn *func;
5723   void            *ctx;
5724 
5725   PetscFunctionBegin;
5726   PetscCall(TSGetRHSJacobian(ts, &J, &B, &func, &ctx));
5727   PetscCall((*func)(ts, 0.0, ts->vec_sol, J, B, ctx));
5728   PetscCall(MatShellTestMult(J, RHSWrapperFunction_TSRHSJacobianTest, ts->vec_sol, ts, flg));
5729   PetscFunctionReturn(PETSC_SUCCESS);
5730 }
5731 
5732 /*@
5733   TSRHSJacobianTestTranspose - Compares the multiply transpose routine provided to the `MATSHELL` with differencing on the `TS` given RHS function.
5734 
5735   Logically Collective
5736 
5737   Input Parameter:
5738 . ts - the time stepping routine
5739 
5740   Output Parameter:
5741 . flg - `PETSC_TRUE` if the multiply is likely correct
5742 
5743   Options Database Key:
5744 . -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view - run the test at each timestep of the integrator
5745 
5746   Level: advanced
5747 
5748   Notes:
5749   This only works for problems defined using `TSSetRHSFunction()` and Jacobian NOT `TSSetIFunction()` and Jacobian
5750 
5751 .seealso: [](ch_ts), `TS`, `Mat`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellTestMultTranspose()`, `TSRHSJacobianTest()`
5752 @*/
5753 PetscErrorCode TSRHSJacobianTestTranspose(TS ts, PetscBool *flg)
5754 {
5755   Mat              J, B;
5756   void            *ctx;
5757   TSRHSJacobianFn *func;
5758 
5759   PetscFunctionBegin;
5760   PetscCall(TSGetRHSJacobian(ts, &J, &B, &func, &ctx));
5761   PetscCall((*func)(ts, 0.0, ts->vec_sol, J, B, ctx));
5762   PetscCall(MatShellTestMultTranspose(J, RHSWrapperFunction_TSRHSJacobianTest, ts->vec_sol, ts, flg));
5763   PetscFunctionReturn(PETSC_SUCCESS);
5764 }
5765 
5766 /*@
5767   TSSetUseSplitRHSFunction - Use the split RHSFunction when a multirate method is used.
5768 
5769   Logically Collective
5770 
5771   Input Parameters:
5772 + ts                   - timestepping context
5773 - use_splitrhsfunction - `PETSC_TRUE` indicates that the split RHSFunction will be used
5774 
5775   Options Database Key:
5776 . -ts_use_splitrhsfunction - <true,false>
5777 
5778   Level: intermediate
5779 
5780   Note:
5781   This is only for multirate methods
5782 
5783 .seealso: [](ch_ts), `TS`, `TSGetUseSplitRHSFunction()`
5784 @*/
5785 PetscErrorCode TSSetUseSplitRHSFunction(TS ts, PetscBool use_splitrhsfunction)
5786 {
5787   PetscFunctionBegin;
5788   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5789   ts->use_splitrhsfunction = use_splitrhsfunction;
5790   PetscFunctionReturn(PETSC_SUCCESS);
5791 }
5792 
5793 /*@
5794   TSGetUseSplitRHSFunction - Gets whether to use the split RHSFunction when a multirate method is used.
5795 
5796   Not Collective
5797 
5798   Input Parameter:
5799 . ts - timestepping context
5800 
5801   Output Parameter:
5802 . use_splitrhsfunction - `PETSC_TRUE` indicates that the split RHSFunction will be used
5803 
5804   Level: intermediate
5805 
5806 .seealso: [](ch_ts), `TS`, `TSSetUseSplitRHSFunction()`
5807 @*/
5808 PetscErrorCode TSGetUseSplitRHSFunction(TS ts, PetscBool *use_splitrhsfunction)
5809 {
5810   PetscFunctionBegin;
5811   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5812   *use_splitrhsfunction = ts->use_splitrhsfunction;
5813   PetscFunctionReturn(PETSC_SUCCESS);
5814 }
5815 
5816 /*@
5817   TSSetMatStructure - sets the relationship between the nonzero structure of the RHS Jacobian matrix to the IJacobian matrix.
5818 
5819   Logically  Collective
5820 
5821   Input Parameters:
5822 + ts  - the time-stepper
5823 - str - the structure (the default is `UNKNOWN_NONZERO_PATTERN`)
5824 
5825   Level: intermediate
5826 
5827   Note:
5828   When the relationship between the nonzero structures is known and supplied the solution process can be much faster
5829 
5830 .seealso: [](ch_ts), `TS`, `MatAXPY()`, `MatStructure`
5831  @*/
5832 PetscErrorCode TSSetMatStructure(TS ts, MatStructure str)
5833 {
5834   PetscFunctionBegin;
5835   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5836   ts->axpy_pattern = str;
5837   PetscFunctionReturn(PETSC_SUCCESS);
5838 }
5839 
5840 /*@
5841   TSSetEvaluationTimes - sets the evaluation points. The solution will be computed and stored for each time requested
5842 
5843   Collective
5844 
5845   Input Parameters:
5846 + ts          - the time-stepper
5847 . n           - number of the time points
5848 - time_points - array of the time points
5849 
5850   Options Database Key:
5851 . -ts_eval_times <t0,...tn> - Sets the evaluation times
5852 
5853   Level: intermediate
5854 
5855   Notes:
5856   The elements in `time_points` must be all increasing. They correspond to the intermediate points for time integration.
5857 
5858   `TS_EXACTFINALTIME_MATCHSTEP` must be used to make the last time step in each sub-interval match the intermediate points specified.
5859 
5860   The intermediate solutions are saved in a vector array that can be accessed with `TSGetEvaluationSolutions()`. Thus using evaluation times may
5861   pressure the memory system when using a large number of time points.
5862 
5863 .seealso: [](ch_ts), `TS`, `TSGetEvaluationTimes()`, `TSGetEvaluationSolutions()`, `TSSetTimeSpan()`
5864  @*/
5865 PetscErrorCode TSSetEvaluationTimes(TS ts, PetscInt n, PetscReal *time_points)
5866 {
5867   PetscBool is_sorted;
5868 
5869   PetscFunctionBegin;
5870   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5871   if (ts->eval_times && n != ts->eval_times->num_time_points) {
5872     PetscCall(PetscFree(ts->eval_times->time_points));
5873     PetscCall(VecDestroyVecs(ts->eval_times->num_time_points, &ts->eval_times->sol_vecs));
5874     PetscCall(PetscMalloc1(n, &ts->eval_times->time_points));
5875   }
5876   if (!ts->eval_times) {
5877     TSEvaluationTimes eval_times;
5878     PetscCall(PetscNew(&eval_times));
5879     PetscCall(PetscMalloc1(n, &eval_times->time_points));
5880     eval_times->reltol  = 1e-6;
5881     eval_times->abstol  = 10 * PETSC_MACHINE_EPSILON;
5882     eval_times->worktol = 0;
5883     ts->eval_times      = eval_times;
5884   }
5885   ts->eval_times->num_time_points = n;
5886   PetscCall(PetscSortedReal(n, time_points, &is_sorted));
5887   PetscCheck(is_sorted, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "time_points array must be sorted");
5888   PetscCall(PetscArraycpy(ts->eval_times->time_points, time_points, n));
5889   PetscFunctionReturn(PETSC_SUCCESS);
5890 }
5891 
5892 /*@C
5893   TSGetEvaluationTimes - gets the evaluation times set with `TSSetEvaluationTimes()`
5894 
5895   Not Collective
5896 
5897   Input Parameter:
5898 . ts - the time-stepper
5899 
5900   Output Parameters:
5901 + n           - number of the time points
5902 - time_points - array of the time points
5903 
5904   Level: beginner
5905 
5906   Note:
5907   The values obtained are valid until the `TS` object is destroyed.
5908 
5909   Both `n` and `time_points` can be `NULL`.
5910 
5911   Also used to see time points set by `TSSetTimeSpan()`.
5912 
5913 .seealso: [](ch_ts), `TS`, `TSSetEvaluationTimes()`, `TSGetEvaluationSolutions()`
5914  @*/
5915 PetscErrorCode TSGetEvaluationTimes(TS ts, PetscInt *n, const PetscReal *time_points[])
5916 {
5917   PetscFunctionBegin;
5918   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5919   if (n) PetscAssertPointer(n, 2);
5920   if (time_points) PetscAssertPointer(time_points, 3);
5921   if (!ts->eval_times) {
5922     if (n) *n = 0;
5923     if (time_points) *time_points = NULL;
5924   } else {
5925     if (n) *n = ts->eval_times->num_time_points;
5926     if (time_points) *time_points = ts->eval_times->time_points;
5927   }
5928   PetscFunctionReturn(PETSC_SUCCESS);
5929 }
5930 
5931 /*@C
5932   TSGetEvaluationSolutions - Get the number of solutions and the solutions at the evaluation time points specified
5933 
5934   Input Parameter:
5935 . ts - the `TS` context obtained from `TSCreate()`
5936 
5937   Output Parameters:
5938 + nsol      - the number of solutions
5939 . sol_times - array of solution times corresponding to the solution vectors. See note below
5940 - Sols      - the solution vectors
5941 
5942   Level: intermediate
5943 
5944   Notes:
5945   Both `nsol` and `Sols` can be `NULL`.
5946 
5947   Some time points in the evaluation points may be skipped by `TS` so that `nsol` is less than the number of points specified by `TSSetEvaluationTimes()`.
5948   For example, manipulating the step size, especially with a reduced precision, may cause `TS` to step over certain evaluation times.
5949 
5950   Also used to see view solutions requested by `TSSetTimeSpan()`.
5951 
5952 .seealso: [](ch_ts), `TS`, `TSSetEvaluationTimes()`, `TSGetEvaluationTimes()`
5953 @*/
5954 PetscErrorCode TSGetEvaluationSolutions(TS ts, PetscInt *nsol, const PetscReal *sol_times[], Vec **Sols)
5955 {
5956   PetscFunctionBegin;
5957   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5958   if (nsol) PetscAssertPointer(nsol, 2);
5959   if (sol_times) PetscAssertPointer(sol_times, 3);
5960   if (Sols) PetscAssertPointer(Sols, 4);
5961   if (!ts->eval_times) {
5962     if (nsol) *nsol = 0;
5963     if (sol_times) *sol_times = NULL;
5964     if (Sols) *Sols = NULL;
5965   } else {
5966     if (nsol) *nsol = ts->eval_times->sol_ctr;
5967     if (sol_times) *sol_times = ts->eval_times->sol_times;
5968     if (Sols) *Sols = ts->eval_times->sol_vecs;
5969   }
5970   PetscFunctionReturn(PETSC_SUCCESS);
5971 }
5972 
5973 /*@
5974   TSSetTimeSpan - sets the time span. The solution will be computed and stored for each time requested in the span
5975 
5976   Collective
5977 
5978   Input Parameters:
5979 + ts         - the time-stepper
5980 . n          - number of the time points (>=2)
5981 - span_times - array of the time points. The first element and the last element are the initial time and the final time respectively.
5982 
5983   Options Database Key:
5984 . -ts_time_span <t0,...tf> - Sets the time span
5985 
5986   Level: intermediate
5987 
5988   Notes:
5989   This function is identical to `TSSetEvaluationTimes()`, except that it also sets the initial time and final time for the `ts` to the first and last `span_times` entries.
5990 
5991   The elements in tspan must be all increasing. They correspond to the intermediate points for time integration.
5992 
5993   `TS_EXACTFINALTIME_MATCHSTEP` must be used to make the last time step in each sub-interval match the intermediate points specified.
5994 
5995   The intermediate solutions are saved in a vector array that can be accessed with `TSGetEvaluationSolutions()`. Thus using time span may
5996   pressure the memory system when using a large number of span points.
5997 
5998 .seealso: [](ch_ts), `TS`, `TSSetEvaluationTimes()`, `TSGetEvaluationTimes()`, `TSGetEvaluationSolutions()`
5999  @*/
6000 PetscErrorCode TSSetTimeSpan(TS ts, PetscInt n, PetscReal *span_times)
6001 {
6002   PetscFunctionBegin;
6003   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
6004   PetscCheck(n >= 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Minimum time span size is 2 but %" PetscInt_FMT " is provided", n);
6005   PetscCall(TSSetEvaluationTimes(ts, n, span_times));
6006   PetscCall(TSSetTime(ts, span_times[0]));
6007   PetscCall(TSSetMaxTime(ts, span_times[n - 1]));
6008   PetscFunctionReturn(PETSC_SUCCESS);
6009 }
6010 
6011 /*@
6012   TSPruneIJacobianColor - Remove nondiagonal zeros in the Jacobian matrix and update the `MatMFFD` coloring information.
6013 
6014   Collective
6015 
6016   Input Parameters:
6017 + ts - the `TS` context
6018 . J  - Jacobian matrix (not altered in this routine)
6019 - B  - newly computed Jacobian matrix to use with preconditioner
6020 
6021   Level: intermediate
6022 
6023   Notes:
6024   This function improves the `MatFDColoring` performance when the Jacobian matrix was over-allocated or contains
6025   many constant zeros entries, which is typically the case when the matrix is generated by a `DM`
6026   and multiple fields are involved.
6027 
6028   Users need to make sure that the Jacobian matrix is properly filled to reflect the sparsity
6029   structure. For `MatFDColoring`, the values of nonzero entries are not important. So one can
6030   usually call `TSComputeIJacobian()` with randomized input vectors to generate a dummy Jacobian.
6031   `TSComputeIJacobian()` should be called before `TSSolve()` but after `TSSetUp()`.
6032 
6033 .seealso: [](ch_ts), `TS`, `MatFDColoring`, `TSComputeIJacobianDefaultColor()`, `MatEliminateZeros()`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`
6034 @*/
6035 PetscErrorCode TSPruneIJacobianColor(TS ts, Mat J, Mat B)
6036 {
6037   MatColoring   mc            = NULL;
6038   ISColoring    iscoloring    = NULL;
6039   MatFDColoring matfdcoloring = NULL;
6040 
6041   PetscFunctionBegin;
6042   /* Generate new coloring after eliminating zeros in the matrix */
6043   PetscCall(MatEliminateZeros(B, PETSC_TRUE));
6044   PetscCall(MatColoringCreate(B, &mc));
6045   PetscCall(MatColoringSetDistance(mc, 2));
6046   PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
6047   PetscCall(MatColoringSetFromOptions(mc));
6048   PetscCall(MatColoringApply(mc, &iscoloring));
6049   PetscCall(MatColoringDestroy(&mc));
6050   /* Replace the old coloring with the new one */
6051   PetscCall(MatFDColoringCreate(B, iscoloring, &matfdcoloring));
6052   PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode (*)(void))SNESTSFormFunction, (void *)ts));
6053   PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
6054   PetscCall(MatFDColoringSetUp(B, iscoloring, matfdcoloring));
6055   PetscCall(PetscObjectCompose((PetscObject)B, "TSMatFDColoring", (PetscObject)matfdcoloring));
6056   PetscCall(PetscObjectDereference((PetscObject)matfdcoloring));
6057   PetscCall(ISColoringDestroy(&iscoloring));
6058   PetscFunctionReturn(PETSC_SUCCESS);
6059 }
6060