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
TSAdaptChoose_History(TSAdapt adapt,TS ts,PetscReal h,PetscInt * next_sc,PetscReal * next_h,PetscBool * accept,PetscReal * wlte,PetscReal * wltea,PetscReal * wlter)8d71ae5a4SJacob 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)
9d71ae5a4SJacob 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;
223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
23d949e4a4SStefano Zampini }
24d949e4a4SStefano Zampini
TSAdaptReset_History(TSAdapt adapt)25d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptReset_History(TSAdapt adapt)
26d71ae5a4SJacob Faibussowitsch {
27d949e4a4SStefano Zampini TSAdapt_History *thadapt = (TSAdapt_History *)adapt->data;
28d949e4a4SStefano Zampini
29d949e4a4SStefano Zampini PetscFunctionBegin;
309566063dSJacob Faibussowitsch PetscCall(TSHistoryDestroy(&thadapt->hist));
313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
32d949e4a4SStefano Zampini }
33d949e4a4SStefano Zampini
TSAdaptDestroy_History(TSAdapt adapt)34d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDestroy_History(TSAdapt adapt)
35d71ae5a4SJacob Faibussowitsch {
36d949e4a4SStefano Zampini PetscFunctionBegin;
379566063dSJacob Faibussowitsch PetscCall(TSAdaptReset_History(adapt));
389566063dSJacob Faibussowitsch PetscCall(PetscFree(adapt->data));
393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
40d949e4a4SStefano Zampini }
41d949e4a4SStefano Zampini
42d949e4a4SStefano Zampini /* this is not public, as TSHistory is not a public object */
TSAdaptHistorySetTSHistory(TSAdapt adapt,TSHistory hist,PetscBool backward)43d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt adapt, TSHistory hist, PetscBool backward)
44d71ae5a4SJacob 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));
533ba16761SJacob Faibussowitsch if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
549566063dSJacob Faibussowitsch PetscCall(TSHistoryGetHistory(hist, &n, (const PetscReal **)&hist_t, NULL, NULL));
559566063dSJacob Faibussowitsch PetscCall(TSAdaptHistorySetHistory(adapt, n, hist_t, backward));
563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
57d949e4a4SStefano Zampini }
58d949e4a4SStefano Zampini
59d949e4a4SStefano Zampini /*@
6075017583SStefano Zampini TSAdaptHistoryGetStep - Gets time and time step for a given step number in the history
6175017583SStefano Zampini
62c3339decSBarry Smith Logically Collective
6375017583SStefano Zampini
6475017583SStefano Zampini Input Parameters:
6575017583SStefano Zampini + adapt - the TSAdapt context
6675017583SStefano Zampini - step - the step number
6775017583SStefano Zampini
6875017583SStefano Zampini Output Parameters:
69*dc9a610eSPierre Jolivet + t - the time corresponding to the requested step (can be `NULL`)
70*dc9a610eSPierre Jolivet - dt - the time step to be taken at the requested step (can be `NULL`)
7175017583SStefano Zampini
7275017583SStefano Zampini Level: advanced
7375017583SStefano Zampini
74bcf0153eSBarry Smith Note:
75bcf0153eSBarry Smith 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()`.
76bcf0153eSBarry Smith
771cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetTrajectory()`, `TSADAPTHISTORY`
7875017583SStefano Zampini @*/
TSAdaptHistoryGetStep(TSAdapt adapt,PetscInt step,PetscReal * t,PetscReal * dt)79d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptHistoryGetStep(TSAdapt adapt, PetscInt step, PetscReal *t, PetscReal *dt)
80d71ae5a4SJacob Faibussowitsch {
8175017583SStefano Zampini TSAdapt_History *thadapt;
8275017583SStefano Zampini PetscBool flg;
8375017583SStefano Zampini
8475017583SStefano Zampini PetscFunctionBegin;
8575017583SStefano Zampini PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
8675017583SStefano Zampini PetscValidLogicalCollectiveInt(adapt, step, 2);
879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
883c633725SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)adapt)->type_name);
8975017583SStefano Zampini thadapt = (TSAdapt_History *)adapt->data;
909566063dSJacob Faibussowitsch PetscCall(TSHistoryGetTimeStep(thadapt->hist, thadapt->bw, step, dt));
919566063dSJacob Faibussowitsch PetscCall(TSHistoryGetTime(thadapt->hist, thadapt->bw, step, t));
923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9375017583SStefano Zampini }
9475017583SStefano Zampini
9575017583SStefano Zampini /*@
96bcf0153eSBarry Smith TSAdaptHistorySetHistory - Sets the time history in the adaptor
97d949e4a4SStefano Zampini
98c3339decSBarry Smith Logically Collective
99d949e4a4SStefano Zampini
100d949e4a4SStefano Zampini Input Parameters:
101bcf0153eSBarry Smith + adapt - the `TSAdapt` context
102d949e4a4SStefano Zampini . n - size of the time history
103d949e4a4SStefano Zampini . hist - the time history
104d949e4a4SStefano Zampini - backward - if the time history has to be followed backward
105d949e4a4SStefano Zampini
106d949e4a4SStefano Zampini Level: advanced
107d949e4a4SStefano Zampini
108bcf0153eSBarry Smith Note:
109bcf0153eSBarry Smith 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()`.
110bcf0153eSBarry Smith
1111cc06b55SBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetTrajectory()`, `TSADAPTHISTORY`
112d949e4a4SStefano Zampini @*/
TSAdaptHistorySetHistory(TSAdapt adapt,PetscInt n,PetscReal hist[],PetscBool backward)113d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptHistorySetHistory(TSAdapt adapt, PetscInt n, PetscReal hist[], PetscBool backward)
114d71ae5a4SJacob Faibussowitsch {
115d949e4a4SStefano Zampini TSAdapt_History *thadapt;
116d949e4a4SStefano Zampini PetscBool flg;
117d949e4a4SStefano Zampini
118d949e4a4SStefano Zampini PetscFunctionBegin;
119d949e4a4SStefano Zampini PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
120d949e4a4SStefano Zampini PetscValidLogicalCollectiveInt(adapt, n, 2);
1214f572ea9SToby Isaac PetscAssertPointer(hist, 3);
122d949e4a4SStefano Zampini PetscValidLogicalCollectiveBool(adapt, backward, 4);
1239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
1243ba16761SJacob Faibussowitsch if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
125d949e4a4SStefano Zampini thadapt = (TSAdapt_History *)adapt->data;
1269566063dSJacob Faibussowitsch PetscCall(TSHistoryDestroy(&thadapt->hist));
1279566063dSJacob Faibussowitsch PetscCall(TSHistoryCreate(PetscObjectComm((PetscObject)adapt), &thadapt->hist));
1289566063dSJacob Faibussowitsch PetscCall(TSHistorySetHistory(thadapt->hist, n, hist, NULL, PETSC_FALSE));
129d949e4a4SStefano Zampini thadapt->bw = backward;
1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
131d949e4a4SStefano Zampini }
132d949e4a4SStefano Zampini
13375017583SStefano Zampini /*@
134bcf0153eSBarry Smith TSAdaptHistorySetTrajectory - Sets a time history in the adaptor from a given `TSTrajectory`
13575017583SStefano Zampini
136c3339decSBarry Smith Logically Collective
13775017583SStefano Zampini
13875017583SStefano Zampini Input Parameters:
139bcf0153eSBarry Smith + adapt - the `TSAdapt` context
140bcf0153eSBarry Smith . tj - the `TSTrajectory` context
14175017583SStefano Zampini - backward - if the time history has to be followed backward
14275017583SStefano Zampini
14375017583SStefano Zampini Level: advanced
14475017583SStefano Zampini
145bcf0153eSBarry Smith Notes:
146bcf0153eSBarry Smith The time history is internally copied, and the user can destroy the `TSTrajectory` if not needed.
147bcf0153eSBarry Smith
148bcf0153eSBarry Smith The user still needs to specify the termination of the solve via `TSSetMaxSteps()`.
149bcf0153eSBarry Smith
1501cc06b55SBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetHistory()`, `TSADAPTHISTORY`, `TSAdapt`
15175017583SStefano Zampini @*/
TSAdaptHistorySetTrajectory(TSAdapt adapt,TSTrajectory tj,PetscBool backward)152d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt adapt, TSTrajectory tj, PetscBool backward)
153d71ae5a4SJacob Faibussowitsch {
15475017583SStefano Zampini PetscBool flg;
15575017583SStefano Zampini
15675017583SStefano Zampini PetscFunctionBegin;
15775017583SStefano Zampini PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
15875017583SStefano Zampini PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 2);
15975017583SStefano Zampini PetscValidLogicalCollectiveBool(adapt, backward, 3);
1609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
1613ba16761SJacob Faibussowitsch if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
1629566063dSJacob Faibussowitsch PetscCall(TSAdaptHistorySetTSHistory(adapt, tj->tsh, backward));
1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16475017583SStefano Zampini }
16575017583SStefano Zampini
166d949e4a4SStefano Zampini /*MC
167d949e4a4SStefano Zampini TSADAPTHISTORY - Time stepping controller that follows a given time history, used for Tangent Linear Model simulations
168d949e4a4SStefano Zampini
169d949e4a4SStefano Zampini Level: developer
170d949e4a4SStefano Zampini
1711cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptHistorySetHistory()`, `TSAdaptType`
172d949e4a4SStefano Zampini M*/
TSAdaptCreate_History(TSAdapt adapt)173d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt adapt)
174d71ae5a4SJacob Faibussowitsch {
175d949e4a4SStefano Zampini TSAdapt_History *thadapt;
176d949e4a4SStefano Zampini
177d949e4a4SStefano Zampini PetscFunctionBegin;
1789566063dSJacob Faibussowitsch PetscCall(PetscNew(&thadapt));
179d949e4a4SStefano Zampini adapt->matchstepfac[0] = PETSC_SMALL; /* prevent from accumulation errors */
180d949e4a4SStefano Zampini adapt->matchstepfac[1] = 0.0; /* we will always match the final step, prevent TSAdaptChoose to mess with it */
181d949e4a4SStefano Zampini adapt->data = thadapt;
182d949e4a4SStefano Zampini
183d949e4a4SStefano Zampini adapt->ops->choose = TSAdaptChoose_History;
184d949e4a4SStefano Zampini adapt->ops->reset = TSAdaptReset_History;
185d949e4a4SStefano Zampini adapt->ops->destroy = TSAdaptDestroy_History;
1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
187d949e4a4SStefano Zampini }
188