xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 7b6bb2c608b6fc6714ef38fda02c2dbb91c82665)
1 /*
2   Code for timestepping with implicit Theta method
3 
4   Notes:
5   This method can be applied to DAE.
6 
7   This method is cast as a 1-stage implicit Runge-Kutta method.
8 
9   Theta | Theta
10   -------------
11         |  1
12 
13   To apply a diagonally implicit RK method to DAE, the stage formula
14 
15   X_i = x + h sum_j a_ij X'_j
16 
17   is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i)
18 */
19 #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
20 
21 typedef struct {
22   Vec       X,Xdot;                   /* Storage for one stage */
23   PetscBool extrapolate;
24   PetscReal Theta;
25   PetscReal shift;
26   PetscReal stage_time;
27 } TS_Theta;
28 
29 #undef __FUNCT__
30 #define __FUNCT__ "TSStep_Theta"
31 static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime)
32 {
33   TS_Theta       *th = (TS_Theta*)ts->data;
34   PetscInt       i,its,lits;
35   PetscErrorCode ierr;
36 
37   PetscFunctionBegin;
38   *steps = -ts->steps;
39   *ptime =  ts->ptime;
40 
41   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
42 
43   for (i=0; i<ts->max_steps; i++) {
44     if (ts->ptime + ts->time_step > ts->max_time) break;
45     ierr = TSPreStep(ts);CHKERRQ(ierr);
46 
47     th->stage_time = ts->ptime + th->Theta*ts->time_step;
48     th->shift = 1./(th->Theta*ts->time_step);
49 
50     if (th->extrapolate) {
51       ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
52     } else {
53       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
54     }
55     ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr);
56     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
57     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
58     ts->nonlinear_its += its; ts->linear_its += lits;
59 
60     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
61     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
62     ts->ptime += ts->time_step;
63     ts->steps++;
64 
65     ierr = TSPostStep(ts);CHKERRQ(ierr);
66     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
67   }
68 
69   *steps += ts->steps;
70   *ptime  = ts->ptime;
71   PetscFunctionReturn(0);
72 }
73 
74 /*------------------------------------------------------------*/
75 #undef __FUNCT__
76 #define __FUNCT__ "TSReset_Theta"
77 static PetscErrorCode TSReset_Theta(TS ts)
78 {
79   TS_Theta       *th = (TS_Theta*)ts->data;
80   PetscErrorCode  ierr;
81 
82   PetscFunctionBegin;
83   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
84   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "TSDestroy_Theta"
90 static PetscErrorCode TSDestroy_Theta(TS ts)
91 {
92   PetscErrorCode  ierr;
93 
94   PetscFunctionBegin;
95   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
96   ierr = PetscFree(ts->data);CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99 
100 /*
101   This defines the nonlinear equation that is to be solved with SNES
102   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
103 */
104 #undef __FUNCT__
105 #define __FUNCT__ "SNESTSFormFunction_Theta"
106 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
107 {
108   TS_Theta       *th = (TS_Theta*)ts->data;
109   PetscErrorCode ierr;
110 
111   PetscFunctionBegin;
112   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr);
113   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr);
114   PetscFunctionReturn(0);
115 }
116 
117 #undef __FUNCT__
118 #define __FUNCT__ "SNESTSFormJacobian_Theta"
119 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
120 {
121   TS_Theta       *th = (TS_Theta*)ts->data;
122   PetscErrorCode ierr;
123 
124   PetscFunctionBegin;
125   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
126   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);CHKERRQ(ierr);
127   PetscFunctionReturn(0);
128 }
129 
130 
131 #undef __FUNCT__
132 #define __FUNCT__ "TSSetUp_Theta"
133 static PetscErrorCode TSSetUp_Theta(TS ts)
134 {
135   TS_Theta       *th = (TS_Theta*)ts->data;
136   PetscErrorCode ierr;
137   Vec            res;
138 
139   PetscFunctionBegin;
140   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
141   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
142   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
143   ierr = VecDuplicate(ts->vec_sol,&res);CHKERRQ(ierr);
144   ierr = SNESSetFunction(ts->snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
145   ierr = VecDestroy(&res);CHKERRQ(ierr);
146   /* This is nasty.  SNESSetFromOptions() is usually called in TSSetFromOptions().  With -snes_mf_operator, it will
147   replace A and we don't want to mess with that.  With -snes_mf, A and B will be replaced as well as the function and
148   context.  Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets
149   the Jacobian user context to snes->funP, it will actually be NULL.  This is not a problem because both snes->funP and
150   snes->jacP should be the TS. */
151   {
152     Mat            A,B;
153     PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
154     void           *ctx;
155     ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr);
156     ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr);
157   }
158   PetscFunctionReturn(0);
159 }
160 /*------------------------------------------------------------*/
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "TSSetFromOptions_Theta"
164 static PetscErrorCode TSSetFromOptions_Theta(TS ts)
165 {
166   TS_Theta       *th = (TS_Theta*)ts->data;
167   PetscErrorCode ierr;
168 
169   PetscFunctionBegin;
170   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
171   {
172     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
173     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
174   }
175   ierr = PetscOptionsTail();CHKERRQ(ierr);
176   PetscFunctionReturn(0);
177 }
178 
179 #undef __FUNCT__
180 #define __FUNCT__ "TSView_Theta"
181 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
182 {
183   TS_Theta       *th = (TS_Theta*)ts->data;
184   PetscBool       iascii;
185   PetscErrorCode  ierr;
186 
187   PetscFunctionBegin;
188   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
189   if (iascii) {
190     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
191     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
192   } else {
193     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name);
194   }
195   PetscFunctionReturn(0);
196 }
197 
198 EXTERN_C_BEGIN
199 #undef __FUNCT__
200 #define __FUNCT__ "TSThetaGetTheta_Theta"
201 PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
202 {
203   TS_Theta *th = (TS_Theta*)ts->data;
204 
205   PetscFunctionBegin;
206   *theta = th->Theta;
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "TSThetaSetTheta_Theta"
212 PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
213 {
214   TS_Theta *th = (TS_Theta*)ts->data;
215 
216   PetscFunctionBegin;
217   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
218   th->Theta = theta;
219   PetscFunctionReturn(0);
220 }
221 EXTERN_C_END
222 
223 /* ------------------------------------------------------------ */
224 /*MC
225       TSTHETA - DAE solver using the implicit Theta method
226 
227   Level: beginner
228 
229 .seealso:  TSCreate(), TS, TSSetType()
230 
231 M*/
232 EXTERN_C_BEGIN
233 #undef __FUNCT__
234 #define __FUNCT__ "TSCreate_Theta"
235 PetscErrorCode  TSCreate_Theta(TS ts)
236 {
237   TS_Theta       *th;
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   ts->ops->reset          = TSReset_Theta;
242   ts->ops->destroy        = TSDestroy_Theta;
243   ts->ops->view           = TSView_Theta;
244   ts->ops->setup          = TSSetUp_Theta;
245   ts->ops->step           = TSStep_Theta;
246   ts->ops->setfromoptions = TSSetFromOptions_Theta;
247   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
248   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
249 
250   ts->problem_type = TS_NONLINEAR;
251   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
252 
253   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
254   ts->data = (void*)th;
255 
256   th->extrapolate = PETSC_FALSE;
257   th->Theta       = 0.5;
258 
259   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
260   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
261   PetscFunctionReturn(0);
262 }
263 EXTERN_C_END
264 
265 #undef __FUNCT__
266 #define __FUNCT__ "TSThetaGetTheta"
267 /*@
268   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
269 
270   Not Collective
271 
272   Input Parameter:
273 .  ts - timestepping context
274 
275   Output Parameter:
276 .  theta - stage abscissa
277 
278   Note:
279   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
280 
281   Level: Advanced
282 
283 .seealso: TSThetaSetTheta()
284 @*/
285 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
286 {
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
291   PetscValidPointer(theta,2);
292   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
293   PetscFunctionReturn(0);
294 }
295 
296 #undef __FUNCT__
297 #define __FUNCT__ "TSThetaSetTheta"
298 /*@
299   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
300 
301   Not Collective
302 
303   Input Parameter:
304 +  ts - timestepping context
305 -  theta - stage abscissa
306 
307   Options Database:
308 .  -ts_theta_theta <theta>
309 
310   Level: Intermediate
311 
312 .seealso: TSThetaGetTheta()
313 @*/
314 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
315 {
316   PetscErrorCode ierr;
317 
318   PetscFunctionBegin;
319   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
320   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
321   PetscFunctionReturn(0);
322 }
323