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