xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 7ba3a57c1a41d0d4ac9568ddc4edec2fa54aa2cd) !
1 /*
2   Code for timestepping with implicit Theta method
3 */
4 #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
5 
6 typedef struct {
7   Vec       X,Xdot;                   /* Storage for one stage */
8   Vec       affine;                   /* Affine vector needed for residual at beginning of step */
9   PetscBool extrapolate;
10   PetscBool endpoint;
11   PetscReal Theta;
12   PetscReal shift;
13   PetscReal stage_time;
14 } TS_Theta;
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "TSStep_Theta"
18 static PetscErrorCode TSStep_Theta(TS ts)
19 {
20   TS_Theta       *th = (TS_Theta*)ts->data;
21   PetscInt       its,lits;
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   ts->time_step = ts->next_time_step;
26   th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
27   th->shift = 1./(th->Theta*ts->time_step);
28 
29   if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
30     ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
31     if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
32     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
33     ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
34   }
35   if (th->extrapolate) {
36     ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
37   } else {
38     ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
39   }
40   ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
41   ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
42   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
43   ts->nonlinear_its += its; ts->linear_its += lits;
44 
45   if (th->endpoint) {
46     ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
47   } else {
48     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
49     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
50   }
51   ts->ptime          += ts->time_step;
52   ts->next_time_step  = ts->time_step;
53   ts->steps++;
54   PetscFunctionReturn(0);
55 }
56 
57 #undef __FUNCT__
58 #define __FUNCT__ "TSInterpolate_Theta"
59 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
60 {
61   TS_Theta       *th = (TS_Theta*)ts->data;
62   PetscReal      alpha = t - ts->ptime;
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
67   if (th->endpoint) alpha *= th->Theta;
68   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
69   PetscFunctionReturn(0);
70 }
71 
72 /*------------------------------------------------------------*/
73 #undef __FUNCT__
74 #define __FUNCT__ "TSReset_Theta"
75 static PetscErrorCode TSReset_Theta(TS ts)
76 {
77   TS_Theta       *th = (TS_Theta*)ts->data;
78   PetscErrorCode  ierr;
79 
80   PetscFunctionBegin;
81   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
82   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
83   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
84   PetscFunctionReturn(0);
85 }
86 
87 #undef __FUNCT__
88 #define __FUNCT__ "TSDestroy_Theta"
89 static PetscErrorCode TSDestroy_Theta(TS ts)
90 {
91   PetscErrorCode  ierr;
92 
93   PetscFunctionBegin;
94   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
95   ierr = PetscFree(ts->data);CHKERRQ(ierr);
96   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
97   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
98   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
99   PetscFunctionReturn(0);
100 }
101 
102 /*
103   This defines the nonlinear equation that is to be solved with SNES
104   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
105 */
106 #undef __FUNCT__
107 #define __FUNCT__ "SNESTSFormFunction_Theta"
108 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
109 {
110   TS_Theta       *th = (TS_Theta*)ts->data;
111   PetscErrorCode ierr;
112 
113   PetscFunctionBegin;
114   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
115   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr);
116   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
117   PetscFunctionReturn(0);
118 }
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "SNESTSFormJacobian_Theta"
122 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
123 {
124   TS_Theta       *th = (TS_Theta*)ts->data;
125   PetscErrorCode ierr;
126 
127   PetscFunctionBegin;
128   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
129   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "TSSetUp_Theta"
136 static PetscErrorCode TSSetUp_Theta(TS ts)
137 {
138   TS_Theta       *th = (TS_Theta*)ts->data;
139   PetscErrorCode ierr;
140 
141   PetscFunctionBegin;
142   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
143   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
144   PetscFunctionReturn(0);
145 }
146 /*------------------------------------------------------------*/
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "TSSetFromOptions_Theta"
150 static PetscErrorCode TSSetFromOptions_Theta(TS ts)
151 {
152   TS_Theta       *th = (TS_Theta*)ts->data;
153   PetscErrorCode ierr;
154 
155   PetscFunctionBegin;
156   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
157   {
158     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
159     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
160     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,PETSC_NULL);CHKERRQ(ierr);
161     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
162   }
163   ierr = PetscOptionsTail();CHKERRQ(ierr);
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "TSView_Theta"
169 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
170 {
171   TS_Theta       *th = (TS_Theta*)ts->data;
172   PetscBool       iascii;
173   PetscErrorCode  ierr;
174 
175   PetscFunctionBegin;
176   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
177   if (iascii) {
178     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
179     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
180   }
181   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
182   PetscFunctionReturn(0);
183 }
184 
185 EXTERN_C_BEGIN
186 #undef __FUNCT__
187 #define __FUNCT__ "TSThetaGetTheta_Theta"
188 PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
189 {
190   TS_Theta *th = (TS_Theta*)ts->data;
191 
192   PetscFunctionBegin;
193   *theta = th->Theta;
194   PetscFunctionReturn(0);
195 }
196 
197 #undef __FUNCT__
198 #define __FUNCT__ "TSThetaSetTheta_Theta"
199 PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
200 {
201   TS_Theta *th = (TS_Theta*)ts->data;
202 
203   PetscFunctionBegin;
204   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
205   th->Theta = theta;
206   PetscFunctionReturn(0);
207 }
208 
209 #undef __FUNCT__
210 #define __FUNCT__ "TSThetaSetEndpoint_Theta"
211 PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
212 {
213   TS_Theta *th = (TS_Theta*)ts->data;
214 
215   PetscFunctionBegin;
216   th->endpoint = flg;
217   PetscFunctionReturn(0);
218 }
219 EXTERN_C_END
220 
221 /* ------------------------------------------------------------ */
222 /*MC
223       TSTHETA - DAE solver using the implicit Theta method
224 
225    Level: beginner
226 
227    Notes:
228    This method can be applied to DAE.
229 
230    This method is cast as a 1-stage implicit Runge-Kutta method.
231 
232 .vb
233   Theta | Theta
234   -------------
235         |  1
236 .ve
237 
238    For the default Theta=0.5, this is also known as the implicit midpoint rule.
239 
240    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
241 
242 .vb
243   0 | 0         0
244   1 | 1-Theta   Theta
245   -------------------
246     | 1-Theta   Theta
247 .ve
248 
249    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
250 
251    To apply a diagonally implicit RK method to DAE, the stage formula
252 
253 $  Y_i = X + h sum_j a_ij Y'_j
254 
255    is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i)
256 
257 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
258 
259 M*/
260 EXTERN_C_BEGIN
261 #undef __FUNCT__
262 #define __FUNCT__ "TSCreate_Theta"
263 PetscErrorCode  TSCreate_Theta(TS ts)
264 {
265   TS_Theta       *th;
266   PetscErrorCode ierr;
267 
268   PetscFunctionBegin;
269   ts->ops->reset          = TSReset_Theta;
270   ts->ops->destroy        = TSDestroy_Theta;
271   ts->ops->view           = TSView_Theta;
272   ts->ops->setup          = TSSetUp_Theta;
273   ts->ops->step           = TSStep_Theta;
274   ts->ops->interpolate    = TSInterpolate_Theta;
275   ts->ops->setfromoptions = TSSetFromOptions_Theta;
276   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
277   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
278 
279   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
280   ts->data = (void*)th;
281 
282   th->extrapolate = PETSC_FALSE;
283   th->Theta       = 0.5;
284 
285   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
286   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
287   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
288   PetscFunctionReturn(0);
289 }
290 EXTERN_C_END
291 
292 #undef __FUNCT__
293 #define __FUNCT__ "TSThetaGetTheta"
294 /*@
295   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
296 
297   Not Collective
298 
299   Input Parameter:
300 .  ts - timestepping context
301 
302   Output Parameter:
303 .  theta - stage abscissa
304 
305   Note:
306   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
307 
308   Level: Advanced
309 
310 .seealso: TSThetaSetTheta()
311 @*/
312 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
313 {
314   PetscErrorCode ierr;
315 
316   PetscFunctionBegin;
317   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
318   PetscValidPointer(theta,2);
319   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
320   PetscFunctionReturn(0);
321 }
322 
323 #undef __FUNCT__
324 #define __FUNCT__ "TSThetaSetTheta"
325 /*@
326   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
327 
328   Not Collective
329 
330   Input Parameter:
331 +  ts - timestepping context
332 -  theta - stage abscissa
333 
334   Options Database:
335 .  -ts_theta_theta <theta>
336 
337   Level: Intermediate
338 
339 .seealso: TSThetaGetTheta()
340 @*/
341 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
342 {
343   PetscErrorCode ierr;
344 
345   PetscFunctionBegin;
346   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
347   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
348   PetscFunctionReturn(0);
349 }
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "TSThetaSetEndpoint"
353 /*@
354   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
355 
356   Not Collective
357 
358   Input Parameter:
359 +  ts - timestepping context
360 -  flg - PETSC_TRUE to use the endpoint variant
361 
362   Options Database:
363 .  -ts_theta_endpoint <flg>
364 
365   Level: Intermediate
366 
367 .seealso: TSTHETA, TSCN
368 @*/
369 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
370 {
371   PetscErrorCode ierr;
372 
373   PetscFunctionBegin;
374   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
375   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
376   PetscFunctionReturn(0);
377 }
378 
379 /*
380  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
381  * The creation functions for these specializations are below.
382  */
383 
384 #undef __FUNCT__
385 #define __FUNCT__ "TSView_BEuler"
386 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
387 {
388   PetscErrorCode ierr;
389 
390   PetscFunctionBegin;
391   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
392   PetscFunctionReturn(0);
393 }
394 
395 /*MC
396       TSBEULER - ODE solver using the implicit backward Euler method
397 
398   Level: beginner
399 
400 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
401 
402 M*/
403 EXTERN_C_BEGIN
404 #undef __FUNCT__
405 #define __FUNCT__ "TSCreate_BEuler"
406 PetscErrorCode  TSCreate_BEuler(TS ts)
407 {
408   PetscErrorCode ierr;
409 
410   PetscFunctionBegin;
411   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
412   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
413   ts->ops->view = TSView_BEuler;
414   PetscFunctionReturn(0);
415 }
416 EXTERN_C_END
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "TSView_CN"
420 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
421 {
422   PetscErrorCode ierr;
423 
424   PetscFunctionBegin;
425   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
426   PetscFunctionReturn(0);
427 }
428 
429 /*MC
430       TSCN - ODE solver using the implicit Crank-Nicolson method.
431 
432   Level: beginner
433 
434   Notes:
435   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
436 
437 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
438 
439 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
440 
441 M*/
442 EXTERN_C_BEGIN
443 #undef __FUNCT__
444 #define __FUNCT__ "TSCreate_CN"
445 PetscErrorCode  TSCreate_CN(TS ts)
446 {
447   PetscErrorCode ierr;
448 
449   PetscFunctionBegin;
450   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
451   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
452   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
453   ts->ops->view = TSView_CN;
454   PetscFunctionReturn(0);
455 }
456 EXTERN_C_END
457