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