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