xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 214bc6a2d8e2ea6a1f725944c2fd032dc3e2c3e7)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown 
4316643e7SJed Brown   Notes:
5316643e7SJed Brown   This method can be applied to DAE.
6316643e7SJed Brown 
7316643e7SJed Brown   This method is cast as a 1-stage implicit Runge-Kutta method.
8316643e7SJed Brown 
9316643e7SJed Brown   Theta | Theta
10316643e7SJed Brown   -------------
11316643e7SJed Brown         |  1
12316643e7SJed Brown 
13316643e7SJed Brown   To apply a diagonally implicit RK method to DAE, the stage formula
14316643e7SJed Brown 
15316643e7SJed Brown   X_i = x + h sum_j a_ij X'_j
16316643e7SJed Brown 
17316643e7SJed Brown   is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i)
18316643e7SJed Brown */
19c6db04a5SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
20316643e7SJed Brown 
21316643e7SJed Brown typedef struct {
22316643e7SJed Brown   Vec       X,Xdot;                   /* Storage for one stage */
23ace3abfcSBarry Smith   PetscBool extrapolate;
24316643e7SJed Brown   PetscReal Theta;
25316643e7SJed Brown   PetscReal shift;
26316643e7SJed Brown   PetscReal stage_time;
27316643e7SJed Brown } TS_Theta;
28316643e7SJed Brown 
29316643e7SJed Brown #undef __FUNCT__
30316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
31316643e7SJed Brown static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime)
32316643e7SJed Brown {
33316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
34186e87acSLisandro Dalcin   PetscInt       i,its,lits;
352b5a38e1SLisandro Dalcin   PetscErrorCode ierr;
36316643e7SJed Brown 
37316643e7SJed Brown   PetscFunctionBegin;
38316643e7SJed Brown   *steps = -ts->steps;
392b5a38e1SLisandro Dalcin   *ptime =  ts->ptime;
402b5a38e1SLisandro Dalcin 
412b5a38e1SLisandro Dalcin   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
42316643e7SJed Brown 
43186e87acSLisandro Dalcin   for (i=0; i<ts->max_steps; i++) {
44316643e7SJed Brown     if (ts->ptime + ts->time_step > ts->max_time) break;
453f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
462b5a38e1SLisandro Dalcin 
47316643e7SJed Brown     th->stage_time = ts->ptime + th->Theta*ts->time_step;
48316643e7SJed Brown     th->shift = 1./(th->Theta*ts->time_step);
49316643e7SJed Brown 
50ace68cafSJed Brown     if (th->extrapolate) {
512b5a38e1SLisandro Dalcin       ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
52ace68cafSJed Brown     } else {
532b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
54ace68cafSJed Brown     }
55316643e7SJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr);
56316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
57316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
58316643e7SJed Brown     ts->nonlinear_its += its; ts->linear_its += lits;
592b5a38e1SLisandro Dalcin 
602b5a38e1SLisandro Dalcin     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
612b5a38e1SLisandro Dalcin     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
622b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
63316643e7SJed Brown     ts->steps++;
642b5a38e1SLisandro Dalcin 
653f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
662b5a38e1SLisandro Dalcin     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
67316643e7SJed Brown   }
68316643e7SJed Brown 
69316643e7SJed Brown   *steps += ts->steps;
70316643e7SJed Brown   *ptime  = ts->ptime;
71316643e7SJed Brown   PetscFunctionReturn(0);
72316643e7SJed Brown }
73316643e7SJed Brown 
74316643e7SJed Brown /*------------------------------------------------------------*/
75316643e7SJed Brown #undef __FUNCT__
76277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
77277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
78316643e7SJed Brown {
79316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
80316643e7SJed Brown   PetscErrorCode  ierr;
81316643e7SJed Brown 
82316643e7SJed Brown   PetscFunctionBegin;
836bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
846bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
85277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
86277b19d0SLisandro Dalcin }
87277b19d0SLisandro Dalcin 
88277b19d0SLisandro Dalcin #undef __FUNCT__
89277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
90277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
91277b19d0SLisandro Dalcin {
92277b19d0SLisandro Dalcin   PetscErrorCode  ierr;
93277b19d0SLisandro Dalcin 
94277b19d0SLisandro Dalcin   PetscFunctionBegin;
95277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
96277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
97316643e7SJed Brown   PetscFunctionReturn(0);
98316643e7SJed Brown }
99316643e7SJed Brown 
100316643e7SJed Brown /*
101316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
1022b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
103316643e7SJed Brown */
104316643e7SJed Brown #undef __FUNCT__
1050f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
1060f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
107316643e7SJed Brown {
108316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
109316643e7SJed Brown   PetscErrorCode ierr;
110316643e7SJed Brown 
111316643e7SJed Brown   PetscFunctionBegin;
1122b5a38e1SLisandro Dalcin   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr);
113*214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
114316643e7SJed Brown   PetscFunctionReturn(0);
115316643e7SJed Brown }
116316643e7SJed Brown 
117316643e7SJed Brown #undef __FUNCT__
1180f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
1190f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
120316643e7SJed Brown {
121316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
122316643e7SJed Brown   PetscErrorCode ierr;
123316643e7SJed Brown 
124316643e7SJed Brown   PetscFunctionBegin;
1250f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
126*214bc6a2SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
127316643e7SJed Brown   PetscFunctionReturn(0);
128316643e7SJed Brown }
129316643e7SJed Brown 
130316643e7SJed Brown 
131316643e7SJed Brown #undef __FUNCT__
132316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
133316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
134316643e7SJed Brown {
135316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
136316643e7SJed Brown   PetscErrorCode ierr;
137316643e7SJed Brown 
138316643e7SJed Brown   PetscFunctionBegin;
139316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
140316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
141316643e7SJed Brown   PetscFunctionReturn(0);
142316643e7SJed Brown }
143316643e7SJed Brown /*------------------------------------------------------------*/
144316643e7SJed Brown 
145316643e7SJed Brown #undef __FUNCT__
146316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
147316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
148316643e7SJed Brown {
149316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
150316643e7SJed Brown   PetscErrorCode ierr;
151316643e7SJed Brown 
152316643e7SJed Brown   PetscFunctionBegin;
153d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
154316643e7SJed Brown   {
155316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
156acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
157316643e7SJed Brown   }
158316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
159316643e7SJed Brown   PetscFunctionReturn(0);
160316643e7SJed Brown }
161316643e7SJed Brown 
162316643e7SJed Brown #undef __FUNCT__
163316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
164316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
165316643e7SJed Brown {
166316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
167ace3abfcSBarry Smith   PetscBool       iascii;
168316643e7SJed Brown   PetscErrorCode  ierr;
169316643e7SJed Brown 
170316643e7SJed Brown   PetscFunctionBegin;
1712692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
172316643e7SJed Brown   if (iascii) {
173316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
174ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
175316643e7SJed Brown   } else {
176e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name);
177316643e7SJed Brown   }
178316643e7SJed Brown   PetscFunctionReturn(0);
179316643e7SJed Brown }
180316643e7SJed Brown 
1810de4c49aSJed Brown EXTERN_C_BEGIN
1820de4c49aSJed Brown #undef __FUNCT__
1830de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
1847087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
1850de4c49aSJed Brown {
1860de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1870de4c49aSJed Brown 
1880de4c49aSJed Brown   PetscFunctionBegin;
1890de4c49aSJed Brown   *theta = th->Theta;
1900de4c49aSJed Brown   PetscFunctionReturn(0);
1910de4c49aSJed Brown }
1920de4c49aSJed Brown 
1930de4c49aSJed Brown #undef __FUNCT__
1940de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
1957087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
1960de4c49aSJed Brown {
1970de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1980de4c49aSJed Brown 
1990de4c49aSJed Brown   PetscFunctionBegin;
200e7be1afaSJed Brown   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
2010de4c49aSJed Brown   th->Theta = theta;
2020de4c49aSJed Brown   PetscFunctionReturn(0);
2030de4c49aSJed Brown }
2040de4c49aSJed Brown EXTERN_C_END
2050de4c49aSJed Brown 
206316643e7SJed Brown /* ------------------------------------------------------------ */
207316643e7SJed Brown /*MC
20896f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
209316643e7SJed Brown 
210316643e7SJed Brown   Level: beginner
211316643e7SJed Brown 
212316643e7SJed Brown .seealso:  TSCreate(), TS, TSSetType()
213316643e7SJed Brown 
214316643e7SJed Brown M*/
215316643e7SJed Brown EXTERN_C_BEGIN
216316643e7SJed Brown #undef __FUNCT__
217316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
2187087cfbeSBarry Smith PetscErrorCode  TSCreate_Theta(TS ts)
219316643e7SJed Brown {
220316643e7SJed Brown   TS_Theta       *th;
221316643e7SJed Brown   PetscErrorCode ierr;
222316643e7SJed Brown 
223316643e7SJed Brown   PetscFunctionBegin;
224277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
225316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
226316643e7SJed Brown   ts->ops->view           = TSView_Theta;
227316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
228316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
229316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
2300f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
2310f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
232316643e7SJed Brown 
233316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
234316643e7SJed Brown   ts->data = (void*)th;
235316643e7SJed Brown 
2366f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
237316643e7SJed Brown   th->Theta       = 0.5;
238316643e7SJed Brown 
2390de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
2400de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
241316643e7SJed Brown   PetscFunctionReturn(0);
242316643e7SJed Brown }
243316643e7SJed Brown EXTERN_C_END
2440de4c49aSJed Brown 
2450de4c49aSJed Brown #undef __FUNCT__
2460de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
2470de4c49aSJed Brown /*@
2480de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
2490de4c49aSJed Brown 
2500de4c49aSJed Brown   Not Collective
2510de4c49aSJed Brown 
2520de4c49aSJed Brown   Input Parameter:
2530de4c49aSJed Brown .  ts - timestepping context
2540de4c49aSJed Brown 
2550de4c49aSJed Brown   Output Parameter:
2560de4c49aSJed Brown .  theta - stage abscissa
2570de4c49aSJed Brown 
2580de4c49aSJed Brown   Note:
2590de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
2600de4c49aSJed Brown 
2610de4c49aSJed Brown   Level: Advanced
2620de4c49aSJed Brown 
2630de4c49aSJed Brown .seealso: TSThetaSetTheta()
2640de4c49aSJed Brown @*/
2657087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
2660de4c49aSJed Brown {
2674ac538c5SBarry Smith   PetscErrorCode ierr;
2680de4c49aSJed Brown 
2690de4c49aSJed Brown   PetscFunctionBegin;
270afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2710de4c49aSJed Brown   PetscValidPointer(theta,2);
2724ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
2730de4c49aSJed Brown   PetscFunctionReturn(0);
2740de4c49aSJed Brown }
2750de4c49aSJed Brown 
2760de4c49aSJed Brown #undef __FUNCT__
2770de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
2780de4c49aSJed Brown /*@
2790de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
2800de4c49aSJed Brown 
2810de4c49aSJed Brown   Not Collective
2820de4c49aSJed Brown 
2830de4c49aSJed Brown   Input Parameter:
2840de4c49aSJed Brown +  ts - timestepping context
2850de4c49aSJed Brown -  theta - stage abscissa
2860de4c49aSJed Brown 
2870de4c49aSJed Brown   Options Database:
2880de4c49aSJed Brown .  -ts_theta_theta <theta>
2890de4c49aSJed Brown 
2900de4c49aSJed Brown   Level: Intermediate
2910de4c49aSJed Brown 
2920de4c49aSJed Brown .seealso: TSThetaGetTheta()
2930de4c49aSJed Brown @*/
2947087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
2950de4c49aSJed Brown {
2964ac538c5SBarry Smith   PetscErrorCode ierr;
2970de4c49aSJed Brown 
2980de4c49aSJed Brown   PetscFunctionBegin;
299afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3004ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
3010de4c49aSJed Brown   PetscFunctionReturn(0);
3020de4c49aSJed Brown }
303