xref: /petsc/src/ts/adapt/impls/history/adapthist.c (revision 373304a1fd22002640631235ee36a21c261ab6d4)
1 #include <petsc/private/tshistoryimpl.h> /*I "petscts.h" I*/
2 
3 typedef struct {
4   TSHistory hist;
5   PetscBool bw;
6 } TSAdapt_History;
7 
TSAdaptChoose_History(TSAdapt adapt,TS ts,PetscReal h,PetscInt * next_sc,PetscReal * next_h,PetscBool * accept,PetscReal * wlte,PetscReal * wltea,PetscReal * wlter)8 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 {
10   PetscInt         step;
11   TSAdapt_History *thadapt = (TSAdapt_History *)adapt->data;
12 
13   PetscFunctionBegin;
14   PetscCheck(thadapt->hist, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ORDER, "Need to call TSAdaptHistorySetHistory() first");
15   PetscCall(TSGetStepNumber(ts, &step));
16   PetscCall(TSHistoryGetTimeStep(thadapt->hist, thadapt->bw, step + 1, next_h));
17   *accept  = PETSC_TRUE;
18   *next_sc = 0;
19   *wlte    = -1;
20   *wltea   = -1;
21   *wlter   = -1;
22   PetscFunctionReturn(PETSC_SUCCESS);
23 }
24 
TSAdaptReset_History(TSAdapt adapt)25 static PetscErrorCode TSAdaptReset_History(TSAdapt adapt)
26 {
27   TSAdapt_History *thadapt = (TSAdapt_History *)adapt->data;
28 
29   PetscFunctionBegin;
30   PetscCall(TSHistoryDestroy(&thadapt->hist));
31   PetscFunctionReturn(PETSC_SUCCESS);
32 }
33 
TSAdaptDestroy_History(TSAdapt adapt)34 static PetscErrorCode TSAdaptDestroy_History(TSAdapt adapt)
35 {
36   PetscFunctionBegin;
37   PetscCall(TSAdaptReset_History(adapt));
38   PetscCall(PetscFree(adapt->data));
39   PetscFunctionReturn(PETSC_SUCCESS);
40 }
41 
42 /* this is not public, as TSHistory is not a public object */
TSAdaptHistorySetTSHistory(TSAdapt adapt,TSHistory hist,PetscBool backward)43 PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt adapt, TSHistory hist, PetscBool backward)
44 {
45   PetscReal *hist_t;
46   PetscInt   n;
47   PetscBool  flg;
48 
49   PetscFunctionBegin;
50   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
51   PetscValidLogicalCollectiveBool(adapt, backward, 3);
52   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
53   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
54   PetscCall(TSHistoryGetHistory(hist, &n, (const PetscReal **)&hist_t, NULL, NULL));
55   PetscCall(TSAdaptHistorySetHistory(adapt, n, hist_t, backward));
56   PetscFunctionReturn(PETSC_SUCCESS);
57 }
58 
59 /*@
60   TSAdaptHistoryGetStep - Gets time and time step for a given step number in the history
61 
62   Logically Collective
63 
64   Input Parameters:
65 + adapt - the TSAdapt context
66 - step  - the step number
67 
68   Output Parameters:
69 + t  - the time corresponding to the requested step (can be `NULL`)
70 - dt - the time step to be taken at the requested step (can be `NULL`)
71 
72   Level: advanced
73 
74   Note:
75   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()`.
76 
77 .seealso: [](ch_ts), `TS`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetTrajectory()`, `TSADAPTHISTORY`
78 @*/
TSAdaptHistoryGetStep(TSAdapt adapt,PetscInt step,PetscReal * t,PetscReal * dt)79 PetscErrorCode TSAdaptHistoryGetStep(TSAdapt adapt, PetscInt step, PetscReal *t, PetscReal *dt)
80 {
81   TSAdapt_History *thadapt;
82   PetscBool        flg;
83 
84   PetscFunctionBegin;
85   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
86   PetscValidLogicalCollectiveInt(adapt, step, 2);
87   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
88   PetscCheck(flg, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)adapt)->type_name);
89   thadapt = (TSAdapt_History *)adapt->data;
90   PetscCall(TSHistoryGetTimeStep(thadapt->hist, thadapt->bw, step, dt));
91   PetscCall(TSHistoryGetTime(thadapt->hist, thadapt->bw, step, t));
92   PetscFunctionReturn(PETSC_SUCCESS);
93 }
94 
95 /*@
96   TSAdaptHistorySetHistory - Sets the time history in the adaptor
97 
98   Logically Collective
99 
100   Input Parameters:
101 + adapt    - the `TSAdapt` context
102 . n        - size of the time history
103 . hist     - the time history
104 - backward - if the time history has to be followed backward
105 
106   Level: advanced
107 
108   Note:
109   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()`.
110 
111 .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetTrajectory()`, `TSADAPTHISTORY`
112 @*/
TSAdaptHistorySetHistory(TSAdapt adapt,PetscInt n,PetscReal hist[],PetscBool backward)113 PetscErrorCode TSAdaptHistorySetHistory(TSAdapt adapt, PetscInt n, PetscReal hist[], PetscBool backward)
114 {
115   TSAdapt_History *thadapt;
116   PetscBool        flg;
117 
118   PetscFunctionBegin;
119   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
120   PetscValidLogicalCollectiveInt(adapt, n, 2);
121   PetscAssertPointer(hist, 3);
122   PetscValidLogicalCollectiveBool(adapt, backward, 4);
123   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
124   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
125   thadapt = (TSAdapt_History *)adapt->data;
126   PetscCall(TSHistoryDestroy(&thadapt->hist));
127   PetscCall(TSHistoryCreate(PetscObjectComm((PetscObject)adapt), &thadapt->hist));
128   PetscCall(TSHistorySetHistory(thadapt->hist, n, hist, NULL, PETSC_FALSE));
129   thadapt->bw = backward;
130   PetscFunctionReturn(PETSC_SUCCESS);
131 }
132 
133 /*@
134   TSAdaptHistorySetTrajectory - Sets a time history in the adaptor from a given `TSTrajectory`
135 
136   Logically Collective
137 
138   Input Parameters:
139 + adapt    - the `TSAdapt` context
140 . tj       - the `TSTrajectory` context
141 - backward - if the time history has to be followed backward
142 
143   Level: advanced
144 
145   Notes:
146   The time history is internally copied, and the user can destroy the `TSTrajectory` if not needed.
147 
148   The user still needs to specify the termination of the solve via `TSSetMaxSteps()`.
149 
150 .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptHistorySetHistory()`, `TSADAPTHISTORY`, `TSAdapt`
151 @*/
TSAdaptHistorySetTrajectory(TSAdapt adapt,TSTrajectory tj,PetscBool backward)152 PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt adapt, TSTrajectory tj, PetscBool backward)
153 {
154   PetscBool flg;
155 
156   PetscFunctionBegin;
157   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
158   PetscValidHeaderSpecific(tj, TSTRAJECTORY_CLASSID, 2);
159   PetscValidLogicalCollectiveBool(adapt, backward, 3);
160   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTHISTORY, &flg));
161   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
162   PetscCall(TSAdaptHistorySetTSHistory(adapt, tj->tsh, backward));
163   PetscFunctionReturn(PETSC_SUCCESS);
164 }
165 
166 /*MC
167    TSADAPTHISTORY - Time stepping controller that follows a given time history, used for Tangent Linear Model simulations
168 
169    Level: developer
170 
171 .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptHistorySetHistory()`, `TSAdaptType`
172 M*/
TSAdaptCreate_History(TSAdapt adapt)173 PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt adapt)
174 {
175   TSAdapt_History *thadapt;
176 
177   PetscFunctionBegin;
178   PetscCall(PetscNew(&thadapt));
179   adapt->matchstepfac[0] = PETSC_SMALL; /* prevent from accumulation errors */
180   adapt->matchstepfac[1] = 0.0;         /* we will always match the final step, prevent TSAdaptChoose to mess with it */
181   adapt->data            = thadapt;
182 
183   adapt->ops->choose  = TSAdaptChoose_History;
184   adapt->ops->reset   = TSAdaptReset_History;
185   adapt->ops->destroy = TSAdaptDestroy_History;
186   PetscFunctionReturn(PETSC_SUCCESS);
187 }
188