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