xref: /petsc/src/ts/trajectory/interface/traj.c (revision 6dd63270497ad23dcf16ae500a87ff2b2a0b7474)
1 #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
2 #include <petsc/private/tshistoryimpl.h>
3 #include <petscdm.h>
4 
5 PetscFunctionList TSTrajectoryList              = NULL;
6 PetscBool         TSTrajectoryRegisterAllCalled = PETSC_FALSE;
7 PetscClassId      TSTRAJECTORY_CLASSID;
8 PetscLogEvent     TSTrajectory_Set, TSTrajectory_Get, TSTrajectory_GetVecs, TSTrajectory_SetUp;
9 
10 /*@C
11   TSTrajectoryRegister - Adds a way of storing trajectories to the `TS` package
12 
13   Not Collective, No Fortran Support
14 
15   Input Parameters:
16 + sname    - the name of a new user-defined creation routine
17 - function - the creation routine itself
18 
19   Level: developer
20 
21   Note:
22   `TSTrajectoryRegister()` may be called multiple times to add several user-defined tses.
23 
24 .seealso: [](ch_ts), `TSTrajectoryRegisterAll()`
25 @*/
26 PetscErrorCode TSTrajectoryRegister(const char sname[], PetscErrorCode (*function)(TSTrajectory, TS))
27 {
28   PetscFunctionBegin;
29   PetscCall(PetscFunctionListAdd(&TSTrajectoryList, sname, function));
30   PetscFunctionReturn(PETSC_SUCCESS);
31 }
32 
33 /*@
34   TSTrajectorySet - Sets a vector of state in the trajectory object
35 
36   Collective
37 
38   Input Parameters:
39 + tj      - the trajectory object
40 . ts      - the time stepper object (optional)
41 . stepnum - the step number
42 . time    - the current time
43 - X       - the current solution
44 
45   Level: developer
46 
47   Note:
48   Usually one does not call this routine, it is called automatically during `TSSolve()`
49 
50 .seealso: [](ch_ts), `TSTrajectorySetUp()`, `TSTrajectoryDestroy()`, `TSTrajectorySetType()`, `TSTrajectorySetVariableNames()`, `TSGetTrajectory()`, `TSTrajectoryGet()`, `TSTrajectoryGetVecs()`
51 @*/
52 PetscErrorCode TSTrajectorySet(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal time, Vec X)
53 {
54   PetscFunctionBegin;
55   if (!tj) PetscFunctionReturn(PETSC_SUCCESS);
56   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
57   if (ts) PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
58   PetscValidLogicalCollectiveInt(tj, stepnum, 3);
59   PetscValidLogicalCollectiveReal(tj, time, 4);
60   PetscValidHeaderSpecific(X, VEC_CLASSID, 5);
61   PetscCheck(tj->setupcalled, PetscObjectComm((PetscObject)tj), PETSC_ERR_ORDER, "TSTrajectorySetUp should be called first");
62   if (tj->monitor) PetscCall(PetscViewerASCIIPrintf(tj->monitor, "TSTrajectorySet: stepnum %" PetscInt_FMT ", time %g (stages %" PetscInt_FMT ")\n", stepnum, (double)time, (PetscInt)!tj->solution_only));
63   PetscCall(PetscLogEventBegin(TSTrajectory_Set, tj, ts, 0, 0));
64   PetscUseTypeMethod(tj, set, ts, stepnum, time, X);
65   PetscCall(PetscLogEventEnd(TSTrajectory_Set, tj, ts, 0, 0));
66   if (tj->usehistory) PetscCall(TSHistoryUpdate(tj->tsh, stepnum, time));
67   if (tj->lag.caching) tj->lag.Udotcached.time = PETSC_MIN_REAL;
68   PetscFunctionReturn(PETSC_SUCCESS);
69 }
70 
71 /*@
72   TSTrajectoryGetNumSteps - Return the number of steps registered in the `TSTrajectory` via `TSTrajectorySet()`.
73 
74   Not Collective.
75 
76   Input Parameter:
77 . tj - the trajectory object
78 
79   Output Parameter:
80 . steps - the number of steps
81 
82   Level: developer
83 
84 .seealso: [](ch_ts), `TS`, `TSTrajectorySet()`
85 @*/
86 PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory tj, PetscInt *steps)
87 {
88   PetscFunctionBegin;
89   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
90   PetscAssertPointer(steps, 2);
91   PetscCall(TSHistoryGetNumSteps(tj->tsh, steps));
92   PetscFunctionReturn(PETSC_SUCCESS);
93 }
94 
95 /*@
96   TSTrajectoryGet - Updates the solution vector of a time stepper object by querying the `TSTrajectory`
97 
98   Collective
99 
100   Input Parameters:
101 + tj      - the trajectory object
102 . ts      - the time stepper object
103 - stepnum - the step number
104 
105   Output Parameter:
106 . time - the time associated with the step number
107 
108   Level: developer
109 
110   Note:
111   Usually one does not call this routine, it is called automatically during `TSSolve()`
112 
113 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSTrajectorySetUp()`, `TSTrajectoryDestroy()`, `TSTrajectorySetType()`, `TSTrajectorySetVariableNames()`, `TSGetTrajectory()`, `TSTrajectorySet()`, `TSTrajectoryGetVecs()`, `TSGetSolution()`
114 @*/
115 PetscErrorCode TSTrajectoryGet(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal *time)
116 {
117   PetscFunctionBegin;
118   PetscCheck(tj, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "TS solver did not save trajectory");
119   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
120   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
121   PetscValidLogicalCollectiveInt(tj, stepnum, 3);
122   PetscAssertPointer(time, 4);
123   PetscCheck(tj->setupcalled, PetscObjectComm((PetscObject)tj), PETSC_ERR_ORDER, "TSTrajectorySetUp should be called first");
124   PetscCheck(stepnum >= 0, PetscObjectComm((PetscObject)tj), PETSC_ERR_PLIB, "Requesting negative step number");
125   if (tj->monitor) {
126     PetscCall(PetscViewerASCIIPrintf(tj->monitor, "TSTrajectoryGet: stepnum %" PetscInt_FMT ", stages %" PetscInt_FMT "\n", stepnum, (PetscInt)!tj->solution_only));
127     PetscCall(PetscViewerFlush(tj->monitor));
128   }
129   PetscCall(PetscLogEventBegin(TSTrajectory_Get, tj, ts, 0, 0));
130   PetscUseTypeMethod(tj, get, ts, stepnum, time);
131   PetscCall(PetscLogEventEnd(TSTrajectory_Get, tj, ts, 0, 0));
132   PetscFunctionReturn(PETSC_SUCCESS);
133 }
134 
135 /*@
136   TSTrajectoryGetVecs - Reconstructs the vector of state and its time derivative using information from the `TSTrajectory` and, possibly, from the `TS`
137 
138   Collective
139 
140   Input Parameters:
141 + tj      - the trajectory object
142 . ts      - the time stepper object (optional)
143 - stepnum - the requested step number
144 
145   Output Parameters:
146 + time - On input time for the step if step number is `PETSC_DECIDE`, on output the time associated with the step number
147 . U    - state vector (can be `NULL`)
148 - Udot - time derivative of state vector (can be `NULL`)
149 
150   Level: developer
151 
152   Notes:
153   If the step number is `PETSC_DECIDE`, the time argument is used to inquire the trajectory.
154   If the requested time does not match any in the trajectory, Lagrangian interpolations are returned.
155 
156 .seealso: [](ch_ts), `TS`, `TSTrajectory`, `TSTrajectorySetUp()`, `TSTrajectoryDestroy()`, `TSTrajectorySetType()`, `TSTrajectorySetVariableNames()`, `TSGetTrajectory()`, `TSTrajectorySet()`, `TSTrajectoryGet()`
157 @*/
158 PetscErrorCode TSTrajectoryGetVecs(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal *time, Vec U, Vec Udot)
159 {
160   PetscFunctionBegin;
161   PetscCheck(tj, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "TS solver did not save trajectory");
162   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
163   if (ts) PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
164   PetscValidLogicalCollectiveInt(tj, stepnum, 3);
165   PetscAssertPointer(time, 4);
166   if (U) PetscValidHeaderSpecific(U, VEC_CLASSID, 5);
167   if (Udot) PetscValidHeaderSpecific(Udot, VEC_CLASSID, 6);
168   if (!U && !Udot) PetscFunctionReturn(PETSC_SUCCESS);
169   PetscCheck(tj->setupcalled, PetscObjectComm((PetscObject)tj), PETSC_ERR_ORDER, "TSTrajectorySetUp should be called first");
170   PetscCall(PetscLogEventBegin(TSTrajectory_GetVecs, tj, ts, 0, 0));
171   if (tj->monitor) {
172     PetscInt pU, pUdot;
173     pU    = U ? 1 : 0;
174     pUdot = Udot ? 1 : 0;
175     PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Requested by GetVecs %" PetscInt_FMT " %" PetscInt_FMT ": stepnum %" PetscInt_FMT ", time %g\n", pU, pUdot, stepnum, (double)*time));
176     PetscCall(PetscViewerFlush(tj->monitor));
177   }
178   if (U && tj->lag.caching) {
179     PetscObjectId    id;
180     PetscObjectState state;
181 
182     PetscCall(PetscObjectStateGet((PetscObject)U, &state));
183     PetscCall(PetscObjectGetId((PetscObject)U, &id));
184     if (stepnum == PETSC_DECIDE) {
185       if (id == tj->lag.Ucached.id && *time == tj->lag.Ucached.time && state == tj->lag.Ucached.state) U = NULL;
186     } else {
187       if (id == tj->lag.Ucached.id && stepnum == tj->lag.Ucached.step && state == tj->lag.Ucached.state) U = NULL;
188     }
189     if (tj->monitor && !U) {
190       PetscCall(PetscViewerASCIIPushTab(tj->monitor));
191       PetscCall(PetscViewerASCIIPrintf(tj->monitor, "State vector cached\n"));
192       PetscCall(PetscViewerASCIIPopTab(tj->monitor));
193       PetscCall(PetscViewerFlush(tj->monitor));
194     }
195   }
196   if (Udot && tj->lag.caching) {
197     PetscObjectId    id;
198     PetscObjectState state;
199 
200     PetscCall(PetscObjectStateGet((PetscObject)Udot, &state));
201     PetscCall(PetscObjectGetId((PetscObject)Udot, &id));
202     if (stepnum == PETSC_DECIDE) {
203       if (id == tj->lag.Udotcached.id && *time == tj->lag.Udotcached.time && state == tj->lag.Udotcached.state) Udot = NULL;
204     } else {
205       if (id == tj->lag.Udotcached.id && stepnum == tj->lag.Udotcached.step && state == tj->lag.Udotcached.state) Udot = NULL;
206     }
207     if (tj->monitor && !Udot) {
208       PetscCall(PetscViewerASCIIPushTab(tj->monitor));
209       PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Derivative vector cached\n"));
210       PetscCall(PetscViewerASCIIPopTab(tj->monitor));
211       PetscCall(PetscViewerFlush(tj->monitor));
212     }
213   }
214   if (!U && !Udot) {
215     PetscCall(PetscLogEventEnd(TSTrajectory_GetVecs, tj, ts, 0, 0));
216     PetscFunctionReturn(PETSC_SUCCESS);
217   }
218 
219   if (stepnum == PETSC_DECIDE || Udot) { /* reverse search for requested time in TSHistory */
220     if (tj->monitor) PetscCall(PetscViewerASCIIPushTab(tj->monitor));
221     /* cached states will be updated in the function */
222     PetscCall(TSTrajectoryReconstruct_Private(tj, ts, *time, U, Udot));
223     if (tj->monitor) {
224       PetscCall(PetscViewerASCIIPopTab(tj->monitor));
225       PetscCall(PetscViewerFlush(tj->monitor));
226     }
227   } else if (U) { /* we were asked to load from stepnum, use TSTrajectoryGet */
228     TS  fakets = ts;
229     Vec U2;
230 
231     /* use a fake TS if ts is missing */
232     if (!ts) {
233       PetscCall(PetscObjectQuery((PetscObject)tj, "__fake_ts", (PetscObject *)&fakets));
234       if (!fakets) {
235         PetscCall(TSCreate(PetscObjectComm((PetscObject)tj), &fakets));
236         PetscCall(PetscObjectCompose((PetscObject)tj, "__fake_ts", (PetscObject)fakets));
237         PetscCall(PetscObjectDereference((PetscObject)fakets));
238         PetscCall(VecDuplicate(U, &U2));
239         PetscCall(TSSetSolution(fakets, U2));
240         PetscCall(PetscObjectDereference((PetscObject)U2));
241       }
242     }
243     PetscCall(TSTrajectoryGet(tj, fakets, stepnum, time));
244     PetscCall(TSGetSolution(fakets, &U2));
245     PetscCall(VecCopy(U2, U));
246     PetscCall(PetscObjectStateGet((PetscObject)U, &tj->lag.Ucached.state));
247     PetscCall(PetscObjectGetId((PetscObject)U, &tj->lag.Ucached.id));
248     tj->lag.Ucached.time = *time;
249     tj->lag.Ucached.step = stepnum;
250   }
251   PetscCall(PetscLogEventEnd(TSTrajectory_GetVecs, tj, ts, 0, 0));
252   PetscFunctionReturn(PETSC_SUCCESS);
253 }
254 
255 /*@
256   TSTrajectoryViewFromOptions - View a `TSTrajectory` based on values in the options database
257 
258   Collective
259 
260   Input Parameters:
261 + A    - the `TSTrajectory` context
262 . obj  - Optional object that provides prefix used for option name
263 - name - command line option
264 
265   Level: intermediate
266 
267 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectoryView`, `PetscObjectViewFromOptions()`, `TSTrajectoryCreate()`
268 @*/
269 PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A, PetscObject obj, const char name[])
270 {
271   PetscFunctionBegin;
272   PetscValidHeaderSpecific(A, TSTRAJECTORY_CLASSID, 1);
273   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
274   PetscFunctionReturn(PETSC_SUCCESS);
275 }
276 
277 /*@
278   TSTrajectoryView - Prints information about the trajectory object
279 
280   Collective
281 
282   Input Parameters:
283 + tj     - the `TSTrajectory` context obtained from `TSTrajectoryCreate()`
284 - viewer - visualization context
285 
286   Options Database Key:
287 . -ts_trajectory_view - calls `TSTrajectoryView()` at end of `TSAdjointStep()`
288 
289   Level: developer
290 
291   Notes:
292   The available visualization contexts include
293 +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
294 -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
295   output where only the first processor opens
296   the file.  All other processors send their
297   data to the first processor to print.
298 
299   The user can open an alternative visualization context with
300   `PetscViewerASCIIOpen()` - output to a specified file.
301 
302 .seealso: [](ch_ts), `TS`, `TSTrajectory`, `PetscViewer`, `PetscViewerASCIIOpen()`
303 @*/
304 PetscErrorCode TSTrajectoryView(TSTrajectory tj, PetscViewer viewer)
305 {
306   PetscBool iascii;
307 
308   PetscFunctionBegin;
309   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
310   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)tj), &viewer));
311   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
312   PetscCheckSameComm(tj, 1, viewer, 2);
313 
314   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
315   if (iascii) {
316     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tj, viewer));
317     PetscCall(PetscViewerASCIIPrintf(viewer, "  total number of recomputations for adjoint calculation = %" PetscInt_FMT "\n", tj->recomps));
318     PetscCall(PetscViewerASCIIPrintf(viewer, "  disk checkpoint reads = %" PetscInt_FMT "\n", tj->diskreads));
319     PetscCall(PetscViewerASCIIPrintf(viewer, "  disk checkpoint writes = %" PetscInt_FMT "\n", tj->diskwrites));
320     PetscCall(PetscViewerASCIIPushTab(viewer));
321     PetscTryTypeMethod(tj, view, viewer);
322     PetscCall(PetscViewerASCIIPopTab(viewer));
323   }
324   PetscFunctionReturn(PETSC_SUCCESS);
325 }
326 
327 /*@C
328   TSTrajectorySetVariableNames - Sets the name of each component in the solution vector so that it may be saved with the trajectory
329 
330   Collective
331 
332   Input Parameters:
333 + ctx   - the trajectory context
334 - names - the names of the components, final string must be `NULL`
335 
336   Level: intermediate
337 
338   Fortran Notes:
339   Fortran interface is not possible because of the string array argument
340 
341 .seealso: [](ch_ts), `TSTrajectory`, `TSGetTrajectory()`
342 @*/
343 PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory ctx, const char *const *names)
344 {
345   PetscFunctionBegin;
346   PetscValidHeaderSpecific(ctx, TSTRAJECTORY_CLASSID, 1);
347   PetscAssertPointer(names, 2);
348   PetscCall(PetscStrArrayDestroy(&ctx->names));
349   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
350   PetscFunctionReturn(PETSC_SUCCESS);
351 }
352 
353 /*@C
354   TSTrajectorySetTransform - Solution vector will be transformed by provided function before being saved to disk
355 
356   Collective
357 
358   Input Parameters:
359 + tj        - the `TSTrajectory` context
360 . transform - the transform function
361 . destroy   - function to destroy the optional context
362 - tctx      - optional context used by transform function
363 
364   Level: intermediate
365 
366 .seealso: [](ch_ts), `TSTrajectorySetVariableNames()`, `TSTrajectory`, `TSMonitorLGSetTransform()`
367 @*/
368 PetscErrorCode TSTrajectorySetTransform(TSTrajectory tj, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
369 {
370   PetscFunctionBegin;
371   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
372   tj->transform        = transform;
373   tj->transformdestroy = destroy;
374   tj->transformctx     = tctx;
375   PetscFunctionReturn(PETSC_SUCCESS);
376 }
377 
378 /*@
379   TSTrajectoryCreate - This function creates an empty trajectory object used to store the time dependent solution of an ODE/DAE
380 
381   Collective
382 
383   Input Parameter:
384 . comm - the communicator
385 
386   Output Parameter:
387 . tj - the trajectory object
388 
389   Level: developer
390 
391   Notes:
392   Usually one does not call this routine, it is called automatically when one calls `TSSetSaveTrajectory()`.
393 
394 .seealso: [](ch_ts), `TS`, `TSTrajectory`, `TSTrajectorySetUp()`, `TSTrajectoryDestroy()`, `TSTrajectorySetType()`, `TSTrajectorySetVariableNames()`, `TSGetTrajectory()`, `TSTrajectorySetKeepFiles()`
395 @*/
396 PetscErrorCode TSTrajectoryCreate(MPI_Comm comm, TSTrajectory *tj)
397 {
398   TSTrajectory t;
399 
400   PetscFunctionBegin;
401   PetscAssertPointer(tj, 2);
402   PetscCall(TSInitializePackage());
403 
404   PetscCall(PetscHeaderCreate(t, TSTRAJECTORY_CLASSID, "TSTrajectory", "Time stepping", "TS", comm, TSTrajectoryDestroy, TSTrajectoryView));
405   t->setupcalled = PETSC_FALSE;
406   PetscCall(TSHistoryCreate(comm, &t->tsh));
407   t->lag.order            = 1;
408   t->lag.L                = NULL;
409   t->lag.T                = NULL;
410   t->lag.W                = NULL;
411   t->lag.WW               = NULL;
412   t->lag.TW               = NULL;
413   t->lag.TT               = NULL;
414   t->lag.caching          = PETSC_TRUE;
415   t->lag.Ucached.id       = 0;
416   t->lag.Ucached.state    = -1;
417   t->lag.Ucached.time     = PETSC_MIN_REAL;
418   t->lag.Ucached.step     = PETSC_INT_MAX;
419   t->lag.Udotcached.id    = 0;
420   t->lag.Udotcached.state = -1;
421   t->lag.Udotcached.time  = PETSC_MIN_REAL;
422   t->lag.Udotcached.step  = PETSC_INT_MAX;
423   t->adjoint_solve_mode   = PETSC_TRUE;
424   t->solution_only        = PETSC_FALSE;
425   t->keepfiles            = PETSC_FALSE;
426   t->usehistory           = PETSC_TRUE;
427   *tj                     = t;
428   PetscCall(TSTrajectorySetFiletemplate(t, "TS-%06" PetscInt_FMT ".bin"));
429   PetscFunctionReturn(PETSC_SUCCESS);
430 }
431 
432 /*@
433   TSTrajectorySetType - Sets the storage method to be used as in a trajectory
434 
435   Collective
436 
437   Input Parameters:
438 + tj   - the `TSTrajectory` context
439 . ts   - the `TS` context
440 - type - a known method
441 
442   Options Database Key:
443 . -ts_trajectory_type <type> - Sets the method; use -help for a list of available methods (for instance, basic)
444 
445   Level: developer
446 
447   Developer Notes:
448   Why does this option require access to the `TS`
449 
450 .seealso: [](ch_ts), `TSTrajectory`, `TS`, `TSTrajectoryCreate()`, `TSTrajectorySetFromOptions()`, `TSTrajectoryDestroy()`, `TSTrajectoryGetType()`
451 @*/
452 PetscErrorCode TSTrajectorySetType(TSTrajectory tj, TS ts, TSTrajectoryType type)
453 {
454   PetscErrorCode (*r)(TSTrajectory, TS);
455   PetscBool match;
456 
457   PetscFunctionBegin;
458   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
459   PetscCall(PetscObjectTypeCompare((PetscObject)tj, type, &match));
460   if (match) PetscFunctionReturn(PETSC_SUCCESS);
461 
462   PetscCall(PetscFunctionListFind(TSTrajectoryList, type, &r));
463   PetscCheck(r, PetscObjectComm((PetscObject)tj), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSTrajectory type: %s", type);
464   PetscTryTypeMethod(tj, destroy);
465   tj->ops->destroy = NULL;
466   PetscCall(PetscMemzero(tj->ops, sizeof(*tj->ops)));
467 
468   PetscCall(PetscObjectChangeTypeName((PetscObject)tj, type));
469   PetscCall((*r)(tj, ts));
470   PetscFunctionReturn(PETSC_SUCCESS);
471 }
472 
473 /*@
474   TSTrajectoryGetType - Gets the trajectory type
475 
476   Collective
477 
478   Input Parameters:
479 + tj - the `TSTrajectory` context
480 - ts - the `TS` context
481 
482   Output Parameter:
483 . type - a known method
484 
485   Level: developer
486 
487 .seealso: [](ch_ts), `TS`, `TSTrajectory`, `TSTrajectoryCreate()`, `TSTrajectorySetFromOptions()`, `TSTrajectoryDestroy()`, `TSTrajectorySetType()`
488 @*/
489 PetscErrorCode TSTrajectoryGetType(TSTrajectory tj, TS ts, TSTrajectoryType *type)
490 {
491   PetscFunctionBegin;
492   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
493   if (type) *type = ((PetscObject)tj)->type_name;
494   PetscFunctionReturn(PETSC_SUCCESS);
495 }
496 
497 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory, TS);
498 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory, TS);
499 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory, TS);
500 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Visualization(TSTrajectory, TS);
501 
502 /*@C
503   TSTrajectoryRegisterAll - Registers all of the `TSTrajectory` storage schecmes in the `TS` package.
504 
505   Not Collective
506 
507   Level: developer
508 
509 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectoryRegister()`
510 @*/
511 PetscErrorCode TSTrajectoryRegisterAll(void)
512 {
513   PetscFunctionBegin;
514   if (TSTrajectoryRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
515   TSTrajectoryRegisterAllCalled = PETSC_TRUE;
516 
517   PetscCall(TSTrajectoryRegister(TSTRAJECTORYBASIC, TSTrajectoryCreate_Basic));
518   PetscCall(TSTrajectoryRegister(TSTRAJECTORYSINGLEFILE, TSTrajectoryCreate_Singlefile));
519   PetscCall(TSTrajectoryRegister(TSTRAJECTORYMEMORY, TSTrajectoryCreate_Memory));
520   PetscCall(TSTrajectoryRegister(TSTRAJECTORYVISUALIZATION, TSTrajectoryCreate_Visualization));
521   PetscFunctionReturn(PETSC_SUCCESS);
522 }
523 
524 /*@
525   TSTrajectoryReset - Resets a trajectory context
526 
527   Collective
528 
529   Input Parameter:
530 . tj - the `TSTrajectory` context obtained from `TSGetTrajectory()`
531 
532   Level: developer
533 
534 .seealso: [](ch_ts), `TS`, `TSTrajectory`, `TSTrajectoryCreate()`, `TSTrajectorySetUp()`
535 @*/
536 PetscErrorCode TSTrajectoryReset(TSTrajectory tj)
537 {
538   PetscFunctionBegin;
539   if (!tj) PetscFunctionReturn(PETSC_SUCCESS);
540   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
541   PetscTryTypeMethod(tj, reset);
542   PetscCall(PetscFree(tj->dirfiletemplate));
543   PetscCall(TSHistoryDestroy(&tj->tsh));
544   PetscCall(TSHistoryCreate(PetscObjectComm((PetscObject)tj), &tj->tsh));
545   tj->setupcalled = PETSC_FALSE;
546   PetscFunctionReturn(PETSC_SUCCESS);
547 }
548 
549 /*@
550   TSTrajectoryDestroy - Destroys a trajectory context
551 
552   Collective
553 
554   Input Parameter:
555 . tj - the `TSTrajectory` context obtained from `TSTrajectoryCreate()`
556 
557   Level: developer
558 
559 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectoryCreate()`, `TSTrajectorySetUp()`
560 @*/
561 PetscErrorCode TSTrajectoryDestroy(TSTrajectory *tj)
562 {
563   PetscFunctionBegin;
564   if (!*tj) PetscFunctionReturn(PETSC_SUCCESS);
565   PetscValidHeaderSpecific(*tj, TSTRAJECTORY_CLASSID, 1);
566   if (--((PetscObject)*tj)->refct > 0) {
567     *tj = NULL;
568     PetscFunctionReturn(PETSC_SUCCESS);
569   }
570 
571   PetscCall(TSTrajectoryReset(*tj));
572   PetscCall(TSHistoryDestroy(&(*tj)->tsh));
573   PetscCall(VecDestroyVecs((*tj)->lag.order + 1, &(*tj)->lag.W));
574   PetscCall(PetscFree5((*tj)->lag.L, (*tj)->lag.T, (*tj)->lag.WW, (*tj)->lag.TT, (*tj)->lag.TW));
575   PetscCall(VecDestroy(&(*tj)->U));
576   PetscCall(VecDestroy(&(*tj)->Udot));
577 
578   if ((*tj)->transformdestroy) PetscCall((*(*tj)->transformdestroy)((*tj)->transformctx));
579   PetscTryTypeMethod(*tj, destroy);
580   if (!((*tj)->keepfiles)) {
581     PetscMPIInt rank;
582     MPI_Comm    comm;
583 
584     PetscCall(PetscObjectGetComm((PetscObject)*tj, &comm));
585     PetscCallMPI(MPI_Comm_rank(comm, &rank));
586     if (rank == 0 && (*tj)->dirname) { /* we own the directory, so we run PetscRMTree on it */
587       PetscCall(PetscRMTree((*tj)->dirname));
588     }
589   }
590   PetscCall(PetscStrArrayDestroy(&(*tj)->names));
591   PetscCall(PetscFree((*tj)->dirname));
592   PetscCall(PetscFree((*tj)->filetemplate));
593   PetscCall(PetscHeaderDestroy(tj));
594   PetscFunctionReturn(PETSC_SUCCESS);
595 }
596 
597 /*
598   TSTrajectorySetTypeFromOptions_Private - Sets the type of `TSTrajectory` from user options.
599 
600   Collective
601 
602   Input Parameter:
603 + tj - the `TSTrajectory` context
604 - ts - the TS context
605 
606   Options Database Key:
607 . -ts_trajectory_type <type> - TSTRAJECTORYBASIC, TSTRAJECTORYMEMORY, TSTRAJECTORYSINGLEFILE, TSTRAJECTORYVISUALIZATION
608 
609   Level: developer
610 
611 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectoryType`, `TSTrajectorySetFromOptions()`, `TSTrajectorySetType()`
612 */
613 static PetscErrorCode TSTrajectorySetTypeFromOptions_Private(PetscOptionItems PetscOptionsObject, TSTrajectory tj, TS ts)
614 {
615   PetscBool   opt;
616   const char *defaultType;
617   char        typeName[256];
618 
619   PetscFunctionBegin;
620   if (((PetscObject)tj)->type_name) defaultType = ((PetscObject)tj)->type_name;
621   else defaultType = TSTRAJECTORYBASIC;
622 
623   PetscCall(TSTrajectoryRegisterAll());
624   PetscCall(PetscOptionsFList("-ts_trajectory_type", "TSTrajectory method", "TSTrajectorySetType", TSTrajectoryList, defaultType, typeName, 256, &opt));
625   if (opt) {
626     PetscCall(TSTrajectorySetType(tj, ts, typeName));
627   } else {
628     PetscCall(TSTrajectorySetType(tj, ts, defaultType));
629   }
630   PetscFunctionReturn(PETSC_SUCCESS);
631 }
632 
633 /*@
634   TSTrajectorySetUseHistory - Use `TSHistory` in `TSTrajectory`
635 
636   Collective
637 
638   Input Parameters:
639 + tj  - the `TSTrajectory` context
640 - flg - `PETSC_TRUE` to save, `PETSC_FALSE` to disable
641 
642   Options Database Key:
643 . -ts_trajectory_use_history - have it use `TSHistory`
644 
645   Level: advanced
646 
647 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`, `TSTrajectorySetUp()`
648 @*/
649 PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory tj, PetscBool flg)
650 {
651   PetscFunctionBegin;
652   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
653   PetscValidLogicalCollectiveBool(tj, flg, 2);
654   tj->usehistory = flg;
655   PetscFunctionReturn(PETSC_SUCCESS);
656 }
657 
658 /*@
659   TSTrajectorySetMonitor - Monitor the schedules generated by the `TSTrajectory` checkpointing controller
660 
661   Collective
662 
663   Input Parameters:
664 + tj  - the `TSTrajectory` context
665 - flg - `PETSC_TRUE` to active a monitor, `PETSC_FALSE` to disable
666 
667   Options Database Key:
668 . -ts_trajectory_monitor - print `TSTrajectory` information
669 
670   Level: developer
671 
672 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`, `TSTrajectorySetUp()`
673 @*/
674 PetscErrorCode TSTrajectorySetMonitor(TSTrajectory tj, PetscBool flg)
675 {
676   PetscFunctionBegin;
677   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
678   PetscValidLogicalCollectiveBool(tj, flg, 2);
679   if (flg) tj->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)tj));
680   else tj->monitor = NULL;
681   PetscFunctionReturn(PETSC_SUCCESS);
682 }
683 
684 /*@
685   TSTrajectorySetKeepFiles - Keep the files generated by the `TSTrajectory` once the program is done
686 
687   Collective
688 
689   Input Parameters:
690 + tj  - the `TSTrajectory` context
691 - flg - `PETSC_TRUE` to save, `PETSC_FALSE` to disable
692 
693   Options Database Key:
694 . -ts_trajectory_keep_files - have it keep the files
695 
696   Level: advanced
697 
698   Note:
699   By default the `TSTrajectory` used for adjoint computations, `TSTRAJECTORYBASIC`, removes the files it generates at the end of the run. This causes the files to be kept.
700 
701 .seealso: [](ch_ts), `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`, `TSTrajectorySetUp()`, `TSTrajectorySetMonitor()`
702 @*/
703 PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory tj, PetscBool flg)
704 {
705   PetscFunctionBegin;
706   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
707   PetscValidLogicalCollectiveBool(tj, flg, 2);
708   tj->keepfiles = flg;
709   PetscFunctionReturn(PETSC_SUCCESS);
710 }
711 
712 /*@
713   TSTrajectorySetDirname - Specify the name of the directory where `TSTrajectory` disk checkpoints are stored.
714 
715   Collective
716 
717   Input Parameters:
718 + tj      - the `TSTrajectory` context
719 - dirname - the directory name
720 
721   Options Database Key:
722 . -ts_trajectory_dirname - set the directory name
723 
724   Level: developer
725 
726   Notes:
727   The final location of the files is determined by dirname/filetemplate where filetemplate was provided by `TSTrajectorySetFiletemplate()`
728 
729   If this is not called `TSTrajectory` selects a unique new name for the directory
730 
731 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectorySetFiletemplate()`, `TSTrajectorySetUp()`
732 @*/
733 PetscErrorCode TSTrajectorySetDirname(TSTrajectory tj, const char dirname[])
734 {
735   PetscBool flg;
736 
737   PetscFunctionBegin;
738   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
739   PetscCall(PetscStrcmp(tj->dirname, dirname, &flg));
740   PetscCheck(flg || !tj->dirfiletemplate, PetscObjectComm((PetscObject)tj), PETSC_ERR_ARG_WRONGSTATE, "Cannot set directoryname after TSTrajectory has been setup");
741   PetscCall(PetscFree(tj->dirname));
742   PetscCall(PetscStrallocpy(dirname, &tj->dirname));
743   PetscFunctionReturn(PETSC_SUCCESS);
744 }
745 
746 /*@
747   TSTrajectorySetFiletemplate - Specify the name template for the files storing `TSTrajectory` checkpoints.
748 
749   Collective
750 
751   Input Parameters:
752 + tj           - the `TSTrajectory` context
753 - filetemplate - the template
754 
755   Options Database Key:
756 . -ts_trajectory_file_template - set the file name template
757 
758   Level: developer
759 
760   Notes:
761   The name template should be of the form, for example filename-%06" PetscInt_FMT ".bin It should not begin with a leading /
762 
763   The final location of the files is determined by dirname/filetemplate where dirname was provided by `TSTrajectorySetDirname()`. The %06" PetscInt_FMT " is replaced by the
764   timestep counter
765 
766 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectorySetDirname()`, `TSTrajectorySetUp()`
767 @*/
768 PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory tj, const char filetemplate[])
769 {
770   const char *ptr = NULL, *ptr2 = NULL;
771 
772   PetscFunctionBegin;
773   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
774   PetscAssertPointer(filetemplate, 2);
775   PetscCheck(!tj->dirfiletemplate, PetscObjectComm((PetscObject)tj), PETSC_ERR_ARG_WRONGSTATE, "Cannot set filetemplate after TSTrajectory has been setup");
776 
777   PetscCheck(filetemplate[0], PetscObjectComm((PetscObject)tj), PETSC_ERR_USER, "-ts_trajectory_file_template requires a file name template, e.g. filename-%%06" PetscInt_FMT ".bin");
778   /* Do some cursory validation of the input. */
779   PetscCall(PetscStrstr(filetemplate, "%", (char **)&ptr));
780   PetscCheck(ptr, PetscObjectComm((PetscObject)tj), PETSC_ERR_USER, "-ts_trajectory_file_template requires a file name template, e.g. filename-%%06" PetscInt_FMT ".bin");
781   for (ptr++; ptr && *ptr; ptr++) {
782     PetscCall(PetscStrchr(PetscInt_FMT "DiouxX", *ptr, (char **)&ptr2));
783     PetscCheck(ptr2 || (*ptr >= '0' && *ptr <= '9'), PetscObjectComm((PetscObject)tj), PETSC_ERR_USER, "Invalid file template argument to -ts_trajectory_file_template, should look like filename-%%06" PetscInt_FMT ".bin");
784     if (ptr2) break;
785   }
786   PetscCall(PetscFree(tj->filetemplate));
787   PetscCall(PetscStrallocpy(filetemplate, &tj->filetemplate));
788   PetscFunctionReturn(PETSC_SUCCESS);
789 }
790 
791 /*@
792   TSTrajectorySetFromOptions - Sets various `TSTrajectory` parameters from user options.
793 
794   Collective
795 
796   Input Parameters:
797 + tj - the `TSTrajectory` context obtained from `TSGetTrajectory()`
798 - ts - the `TS` context
799 
800   Options Database Keys:
801 + -ts_trajectory_type <type>             - basic, memory, singlefile, visualization
802 . -ts_trajectory_keep_files <true,false> - keep the files generated by the code after the program ends. This is true by default for singlefile and visualization
803 - -ts_trajectory_monitor                 - print `TSTrajectory` information
804 
805   Level: developer
806 
807   Note:
808   This is not normally called directly by users
809 
810 .seealso: [](ch_ts), `TSTrajectory`, `TSSetSaveTrajectory()`, `TSTrajectorySetUp()`
811 @*/
812 PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory tj, TS ts)
813 {
814   PetscBool set, flg;
815   char      dirname[PETSC_MAX_PATH_LEN], filetemplate[PETSC_MAX_PATH_LEN];
816 
817   PetscFunctionBegin;
818   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
819   if (ts) PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
820   PetscObjectOptionsBegin((PetscObject)tj);
821   PetscCall(TSTrajectorySetTypeFromOptions_Private(PetscOptionsObject, tj, ts));
822   PetscCall(PetscOptionsBool("-ts_trajectory_use_history", "Turn on/off usage of TSHistory", NULL, tj->usehistory, &tj->usehistory, NULL));
823   PetscCall(PetscOptionsBool("-ts_trajectory_monitor", "Print checkpointing schedules", "TSTrajectorySetMonitor", tj->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set));
824   if (set) PetscCall(TSTrajectorySetMonitor(tj, flg));
825   PetscCall(PetscOptionsInt("-ts_trajectory_reconstruction_order", "Interpolation order for reconstruction", NULL, tj->lag.order, &tj->lag.order, NULL));
826   PetscCall(PetscOptionsBool("-ts_trajectory_reconstruction_caching", "Turn on/off caching of TSTrajectoryGetVecs input", NULL, tj->lag.caching, &tj->lag.caching, NULL));
827   PetscCall(PetscOptionsBool("-ts_trajectory_adjointmode", "Instruct the trajectory that will be used in a TSAdjointSolve()", NULL, tj->adjoint_solve_mode, &tj->adjoint_solve_mode, NULL));
828   PetscCall(PetscOptionsBool("-ts_trajectory_solution_only", "Checkpoint solution only", "TSTrajectorySetSolutionOnly", tj->solution_only, &tj->solution_only, NULL));
829   PetscCall(PetscOptionsBool("-ts_trajectory_keep_files", "Keep any trajectory files generated during the run", "TSTrajectorySetKeepFiles", tj->keepfiles, &flg, &set));
830   if (set) PetscCall(TSTrajectorySetKeepFiles(tj, flg));
831 
832   PetscCall(PetscOptionsString("-ts_trajectory_dirname", "Directory name for TSTrajectory file", "TSTrajectorySetDirname", NULL, dirname, sizeof(dirname) - 14, &set));
833   if (set) PetscCall(TSTrajectorySetDirname(tj, dirname));
834 
835   PetscCall(PetscOptionsString("-ts_trajectory_file_template", "Template for TSTrajectory file name, use filename-%06" PetscInt_FMT ".bin", "TSTrajectorySetFiletemplate", NULL, filetemplate, sizeof(filetemplate), &set));
836   if (set) PetscCall(TSTrajectorySetFiletemplate(tj, filetemplate));
837 
838   /* Handle specific TSTrajectory options */
839   PetscTryTypeMethod(tj, setfromoptions, PetscOptionsObject);
840   PetscOptionsEnd();
841   PetscFunctionReturn(PETSC_SUCCESS);
842 }
843 
844 /*@
845   TSTrajectorySetUp - Sets up the internal data structures, e.g. stacks, for the later use
846   of a `TS` `TSTrajectory`.
847 
848   Collective
849 
850   Input Parameters:
851 + tj - the `TSTrajectory` context
852 - ts - the TS context obtained from `TSCreate()`
853 
854   Level: developer
855 
856 .seealso: [](ch_ts), `TSTrajectory`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
857 @*/
858 PetscErrorCode TSTrajectorySetUp(TSTrajectory tj, TS ts)
859 {
860   size_t s1, s2;
861 
862   PetscFunctionBegin;
863   if (!tj) PetscFunctionReturn(PETSC_SUCCESS);
864   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
865   if (ts) PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
866   if (tj->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
867 
868   PetscCall(PetscLogEventBegin(TSTrajectory_SetUp, tj, ts, 0, 0));
869   if (!((PetscObject)tj)->type_name) PetscCall(TSTrajectorySetType(tj, ts, TSTRAJECTORYBASIC));
870   PetscTryTypeMethod(tj, setup, ts);
871 
872   tj->setupcalled = PETSC_TRUE;
873 
874   /* Set the counters to zero */
875   tj->recomps    = 0;
876   tj->diskreads  = 0;
877   tj->diskwrites = 0;
878   PetscCall(PetscStrlen(tj->dirname, &s1));
879   PetscCall(PetscStrlen(tj->filetemplate, &s2));
880   PetscCall(PetscFree(tj->dirfiletemplate));
881   PetscCall(PetscMalloc((s1 + s2 + 10) * sizeof(char), &tj->dirfiletemplate));
882   PetscCall(PetscSNPrintf(tj->dirfiletemplate, s1 + s2 + 10, "%s/%s", tj->dirname, tj->filetemplate));
883   PetscCall(PetscLogEventEnd(TSTrajectory_SetUp, tj, ts, 0, 0));
884   PetscFunctionReturn(PETSC_SUCCESS);
885 }
886 
887 /*@
888   TSTrajectorySetSolutionOnly - Tells the trajectory to store just the solution, and not any intermediate stage information
889 
890   Collective
891 
892   Input Parameters:
893 + tj            - the `TSTrajectory` context obtained with `TSGetTrajectory()`
894 - solution_only - the boolean flag
895 
896   Level: developer
897 
898 .seealso: [](ch_ts), `TSTrajectory`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`, `TSTrajectoryGetSolutionOnly()`
899 @*/
900 PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj, PetscBool solution_only)
901 {
902   PetscFunctionBegin;
903   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
904   PetscValidLogicalCollectiveBool(tj, solution_only, 2);
905   tj->solution_only = solution_only;
906   PetscFunctionReturn(PETSC_SUCCESS);
907 }
908 
909 /*@
910   TSTrajectoryGetSolutionOnly - Gets the value set with `TSTrajectorySetSolutionOnly()`.
911 
912   Logically Collective
913 
914   Input Parameter:
915 . tj - the `TSTrajectory` context
916 
917   Output Parameter:
918 . solution_only - the boolean flag
919 
920   Level: developer
921 
922 .seealso: [](ch_ts), `TSTrajectory`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`, `TSTrajectorySetSolutionOnly()`
923 @*/
924 PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory tj, PetscBool *solution_only)
925 {
926   PetscFunctionBegin;
927   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
928   PetscAssertPointer(solution_only, 2);
929   *solution_only = tj->solution_only;
930   PetscFunctionReturn(PETSC_SUCCESS);
931 }
932 
933 /*@
934   TSTrajectoryGetUpdatedHistoryVecs - Get updated state and time-derivative history vectors.
935 
936   Collective
937 
938   Input Parameters:
939 + tj   - the `TSTrajectory` context
940 . ts   - the `TS` solver context
941 - time - the requested time
942 
943   Output Parameters:
944 + U    - state vector at given time (can be interpolated)
945 - Udot - time-derivative vector at given time (can be interpolated)
946 
947   Level: developer
948 
949   Notes:
950   The vectors are interpolated if time does not match any time step stored in the `TSTrajectory()`. Pass NULL to not request a vector.
951 
952   This function differs from `TSTrajectoryGetVecs()` since the vectors obtained cannot be modified, and they need to be returned by
953   calling `TSTrajectoryRestoreUpdatedHistoryVecs()`.
954 
955 .seealso: [](ch_ts), `TSTrajectory`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`, `TSTrajectoryRestoreUpdatedHistoryVecs()`, `TSTrajectoryGetVecs()`
956 @*/
957 PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory tj, TS ts, PetscReal time, Vec *U, Vec *Udot)
958 {
959   PetscFunctionBegin;
960   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
961   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
962   PetscValidLogicalCollectiveReal(tj, time, 3);
963   if (U) PetscAssertPointer(U, 4);
964   if (Udot) PetscAssertPointer(Udot, 5);
965   if (U && !tj->U) {
966     DM dm;
967 
968     PetscCall(TSGetDM(ts, &dm));
969     PetscCall(DMCreateGlobalVector(dm, &tj->U));
970   }
971   if (Udot && !tj->Udot) {
972     DM dm;
973 
974     PetscCall(TSGetDM(ts, &dm));
975     PetscCall(DMCreateGlobalVector(dm, &tj->Udot));
976   }
977   PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &time, U ? tj->U : NULL, Udot ? tj->Udot : NULL));
978   if (U) {
979     PetscCall(VecLockReadPush(tj->U));
980     *U = tj->U;
981   }
982   if (Udot) {
983     PetscCall(VecLockReadPush(tj->Udot));
984     *Udot = tj->Udot;
985   }
986   PetscFunctionReturn(PETSC_SUCCESS);
987 }
988 
989 /*@
990   TSTrajectoryRestoreUpdatedHistoryVecs - Restores updated state and time-derivative history vectors obtained with `TSTrajectoryGetUpdatedHistoryVecs()`.
991 
992   Collective
993 
994   Input Parameters:
995 + tj   - the `TSTrajectory` context
996 . U    - state vector at given time (can be interpolated)
997 - Udot - time-derivative vector at given time (can be interpolated)
998 
999   Level: developer
1000 
1001 .seealso: [](ch_ts), `TSTrajectory`, `TSTrajectoryGetUpdatedHistoryVecs()`
1002 @*/
1003 PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory tj, Vec *U, Vec *Udot)
1004 {
1005   PetscFunctionBegin;
1006   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 1);
1007   if (U) PetscValidHeaderSpecific(*U, VEC_CLASSID, 2);
1008   if (Udot) PetscValidHeaderSpecific(*Udot, VEC_CLASSID, 3);
1009   PetscCheck(!U || *U == tj->U, PetscObjectComm((PetscObject)*U), PETSC_ERR_USER, "U was not obtained from TSTrajectoryGetUpdatedHistoryVecs()");
1010   PetscCheck(!Udot || *Udot == tj->Udot, PetscObjectComm((PetscObject)*Udot), PETSC_ERR_USER, "Udot was not obtained from TSTrajectoryGetUpdatedHistoryVecs()");
1011   if (U) {
1012     PetscCall(VecLockReadPop(tj->U));
1013     *U = NULL;
1014   }
1015   if (Udot) {
1016     PetscCall(VecLockReadPop(tj->Udot));
1017     *Udot = NULL;
1018   }
1019   PetscFunctionReturn(PETSC_SUCCESS);
1020 }
1021