xref: /petsc/src/ts/adapt/impls/history/adapthist.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1d949e4a4SStefano Zampini #include <petsc/private/tshistoryimpl.h> /*I "petscts.h" I*/
2d949e4a4SStefano Zampini 
3d949e4a4SStefano Zampini typedef struct {
4d949e4a4SStefano Zampini   TSHistory hist;
5d949e4a4SStefano Zampini   PetscBool bw;
6d949e4a4SStefano Zampini } TSAdapt_History;
7d949e4a4SStefano Zampini 
8*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptChoose_History(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
9*d71ae5a4SJacob Faibussowitsch {
10d949e4a4SStefano Zampini   PetscInt         step;
11d949e4a4SStefano Zampini   TSAdapt_History *thadapt = (TSAdapt_History *)adapt->data;
12d949e4a4SStefano Zampini 
13d949e4a4SStefano Zampini   PetscFunctionBegin;
143c633725SBarry Smith   PetscCheck(thadapt->hist, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ORDER, "Need to call TSAdaptHistorySetHistory() first");
159566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &step));
169566063dSJacob Faibussowitsch   PetscCall(TSHistoryGetTimeStep(thadapt->hist, thadapt->bw, step + 1, next_h));
17d949e4a4SStefano Zampini   *accept  = PETSC_TRUE;
18d949e4a4SStefano Zampini   *next_sc = 0;
19d949e4a4SStefano Zampini   *wlte    = -1;
20d949e4a4SStefano Zampini   *wltea   = -1;
21d949e4a4SStefano Zampini   *wlter   = -1;
22d949e4a4SStefano Zampini   PetscFunctionReturn(0);
23d949e4a4SStefano Zampini }
24d949e4a4SStefano Zampini 
25*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptReset_History(TSAdapt adapt)
26*d71ae5a4SJacob Faibussowitsch {
27d949e4a4SStefano Zampini   TSAdapt_History *thadapt = (TSAdapt_History *)adapt->data;
28d949e4a4SStefano Zampini 
29d949e4a4SStefano Zampini   PetscFunctionBegin;
309566063dSJacob Faibussowitsch   PetscCall(TSHistoryDestroy(&thadapt->hist));
31d949e4a4SStefano Zampini   PetscFunctionReturn(0);
32d949e4a4SStefano Zampini }
33d949e4a4SStefano Zampini 
34*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDestroy_History(TSAdapt adapt)
35*d71ae5a4SJacob Faibussowitsch {
36d949e4a4SStefano Zampini   PetscFunctionBegin;
379566063dSJacob Faibussowitsch   PetscCall(TSAdaptReset_History(adapt));
389566063dSJacob Faibussowitsch   PetscCall(PetscFree(adapt->data));
39d949e4a4SStefano Zampini   PetscFunctionReturn(0);
40d949e4a4SStefano Zampini }
41d949e4a4SStefano Zampini 
42d949e4a4SStefano Zampini /* this is not public, as TSHistory is not a public object */
43*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt adapt, TSHistory hist, PetscBool backward)
44*d71ae5a4SJacob Faibussowitsch {
45d949e4a4SStefano Zampini   PetscReal *hist_t;
46d949e4a4SStefano Zampini   PetscInt   n;
47d949e4a4SStefano Zampini   PetscBool  flg;
48d949e4a4SStefano Zampini 
49d949e4a4SStefano Zampini   PetscFunctionBegin;
50d949e4a4SStefano Zampini   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
51d949e4a4SStefano Zampini   PetscValidLogicalCollectiveBool(adapt, backward, 3);
529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
53d949e4a4SStefano Zampini   if (!flg) PetscFunctionReturn(0);
549566063dSJacob Faibussowitsch   PetscCall(TSHistoryGetHistory(hist, &n, (const PetscReal **)&hist_t, NULL, NULL));
559566063dSJacob Faibussowitsch   PetscCall(TSAdaptHistorySetHistory(adapt, n, hist_t, backward));
56d949e4a4SStefano Zampini   PetscFunctionReturn(0);
57d949e4a4SStefano Zampini }
58d949e4a4SStefano Zampini 
59d949e4a4SStefano Zampini /*@
6075017583SStefano Zampini    TSAdaptHistoryGetStep - Gets time and time step for a given step number in the history
6175017583SStefano Zampini 
6275017583SStefano Zampini    Logically Collective on TSAdapt
6375017583SStefano Zampini 
6475017583SStefano Zampini    Input Parameters:
6575017583SStefano Zampini +  adapt    - the TSAdapt context
6675017583SStefano Zampini -  step     - the step number
6775017583SStefano Zampini 
6875017583SStefano Zampini    Output Parameters:
6975017583SStefano Zampini +  t  - the time corresponding to the requested step (can be NULL)
7075017583SStefano Zampini -  dt - the time step to be taken at the requested step (can be NULL)
7175017583SStefano Zampini 
7275017583SStefano Zampini    Notes: The time history is internally copied, and the user can free the hist array. The user still needs to specify the termination of the solve via TSSetMaxSteps().
7375017583SStefano Zampini 
7475017583SStefano Zampini    Level: advanced
7575017583SStefano Zampini 
76db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetTrajectory()`, `TSADAPTHISTORY`
7775017583SStefano Zampini @*/
78*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptHistoryGetStep(TSAdapt adapt, PetscInt step, PetscReal *t, PetscReal *dt)
79*d71ae5a4SJacob Faibussowitsch {
8075017583SStefano Zampini   TSAdapt_History *thadapt;
8175017583SStefano Zampini   PetscBool        flg;
8275017583SStefano Zampini 
8375017583SStefano Zampini   PetscFunctionBegin;
8475017583SStefano Zampini   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
8575017583SStefano Zampini   PetscValidLogicalCollectiveInt(adapt, step, 2);
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
873c633725SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)adapt)->type_name);
8875017583SStefano Zampini   thadapt = (TSAdapt_History *)adapt->data;
899566063dSJacob Faibussowitsch   PetscCall(TSHistoryGetTimeStep(thadapt->hist, thadapt->bw, step, dt));
909566063dSJacob Faibussowitsch   PetscCall(TSHistoryGetTime(thadapt->hist, thadapt->bw, step, t));
9175017583SStefano Zampini   PetscFunctionReturn(0);
9275017583SStefano Zampini }
9375017583SStefano Zampini 
9475017583SStefano Zampini /*@
95d949e4a4SStefano Zampini    TSAdaptHistorySetHistory - Sets a time history in the adaptor
96d949e4a4SStefano Zampini 
97d949e4a4SStefano Zampini    Logically Collective on TSAdapt
98d949e4a4SStefano Zampini 
99d949e4a4SStefano Zampini    Input Parameters:
100d949e4a4SStefano Zampini +  adapt    - the TSAdapt context
101d949e4a4SStefano Zampini .  n        - size of the time history
102d949e4a4SStefano Zampini .  hist     - the time history
103d949e4a4SStefano Zampini -  backward - if the time history has to be followed backward
104d949e4a4SStefano Zampini 
105d949e4a4SStefano Zampini    Notes: The time history is internally copied, and the user can free the hist array. The user still needs to specify the termination of the solve via TSSetMaxSteps().
106d949e4a4SStefano Zampini 
107d949e4a4SStefano Zampini    Level: advanced
108d949e4a4SStefano Zampini 
109db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetTrajectory()`, `TSADAPTHISTORY`
110d949e4a4SStefano Zampini @*/
111*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptHistorySetHistory(TSAdapt adapt, PetscInt n, PetscReal hist[], PetscBool backward)
112*d71ae5a4SJacob Faibussowitsch {
113d949e4a4SStefano Zampini   TSAdapt_History *thadapt;
114d949e4a4SStefano Zampini   PetscBool        flg;
115d949e4a4SStefano Zampini 
116d949e4a4SStefano Zampini   PetscFunctionBegin;
117d949e4a4SStefano Zampini   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
118d949e4a4SStefano Zampini   PetscValidLogicalCollectiveInt(adapt, n, 2);
119d949e4a4SStefano Zampini   PetscValidRealPointer(hist, 3);
120d949e4a4SStefano Zampini   PetscValidLogicalCollectiveBool(adapt, backward, 4);
1219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
122d949e4a4SStefano Zampini   if (!flg) PetscFunctionReturn(0);
123d949e4a4SStefano Zampini   thadapt = (TSAdapt_History *)adapt->data;
1249566063dSJacob Faibussowitsch   PetscCall(TSHistoryDestroy(&thadapt->hist));
1259566063dSJacob Faibussowitsch   PetscCall(TSHistoryCreate(PetscObjectComm((PetscObject)adapt), &thadapt->hist));
1269566063dSJacob Faibussowitsch   PetscCall(TSHistorySetHistory(thadapt->hist, n, hist, NULL, PETSC_FALSE));
127d949e4a4SStefano Zampini   thadapt->bw = backward;
128d949e4a4SStefano Zampini   PetscFunctionReturn(0);
129d949e4a4SStefano Zampini }
130d949e4a4SStefano Zampini 
13175017583SStefano Zampini /*@
13275017583SStefano Zampini    TSAdaptHistorySetTrajectory - Sets a time history in the adaptor from a given TSTrajectory
13375017583SStefano Zampini 
13475017583SStefano Zampini    Logically Collective on TSAdapt
13575017583SStefano Zampini 
13675017583SStefano Zampini    Input Parameters:
13775017583SStefano Zampini +  adapt    - the TSAdapt context
13875017583SStefano Zampini .  tj       - the TSTrajectory context
13975017583SStefano Zampini -  backward - if the time history has to be followed backward
14075017583SStefano Zampini 
14175017583SStefano Zampini    Notes: The time history is internally copied, and the user can destroy the TSTrajectory if not needed. The user still needs to specify the termination of the solve via TSSetMaxSteps().
14275017583SStefano Zampini 
14375017583SStefano Zampini    Level: advanced
14475017583SStefano Zampini 
145db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetHistory()`, `TSADAPTHISTORY`
14675017583SStefano Zampini @*/
147*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt adapt, TSTrajectory tj, PetscBool backward)
148*d71ae5a4SJacob Faibussowitsch {
14975017583SStefano Zampini   PetscBool flg;
15075017583SStefano Zampini 
15175017583SStefano Zampini   PetscFunctionBegin;
15275017583SStefano Zampini   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
15375017583SStefano Zampini   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 2);
15475017583SStefano Zampini   PetscValidLogicalCollectiveBool(adapt, backward, 3);
1559566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
15675017583SStefano Zampini   if (!flg) PetscFunctionReturn(0);
1579566063dSJacob Faibussowitsch   PetscCall(TSAdaptHistorySetTSHistory(adapt, tj->tsh, backward));
15875017583SStefano Zampini   PetscFunctionReturn(0);
15975017583SStefano Zampini }
16075017583SStefano Zampini 
161d949e4a4SStefano Zampini /*MC
162d949e4a4SStefano Zampini    TSADAPTHISTORY - Time stepping controller that follows a given time history, used for Tangent Linear Model simulations
163d949e4a4SStefano Zampini 
164d949e4a4SStefano Zampini    Level: developer
165d949e4a4SStefano Zampini 
166db781477SPatrick Sanan .seealso: `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptHistorySetHistory()`
167d949e4a4SStefano Zampini M*/
168*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt adapt)
169*d71ae5a4SJacob Faibussowitsch {
170d949e4a4SStefano Zampini   TSAdapt_History *thadapt;
171d949e4a4SStefano Zampini 
172d949e4a4SStefano Zampini   PetscFunctionBegin;
1739566063dSJacob Faibussowitsch   PetscCall(PetscNew(&thadapt));
174d949e4a4SStefano Zampini   adapt->matchstepfac[0] = PETSC_SMALL; /* prevent from accumulation errors */
175d949e4a4SStefano Zampini   adapt->matchstepfac[1] = 0.0;         /* we will always match the final step, prevent TSAdaptChoose to mess with it */
176d949e4a4SStefano Zampini   adapt->data            = thadapt;
177d949e4a4SStefano Zampini 
178d949e4a4SStefano Zampini   adapt->ops->choose  = TSAdaptChoose_History;
179d949e4a4SStefano Zampini   adapt->ops->reset   = TSAdaptReset_History;
180d949e4a4SStefano Zampini   adapt->ops->destroy = TSAdaptDestroy_History;
181d949e4a4SStefano Zampini   PetscFunctionReturn(0);
182d949e4a4SStefano Zampini }
183