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