xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 0f5c6efe1f3f3c860ea310cfc5c66eb104cb8809)
1316643e7SJed Brown #define PETSCTS_DLL
2316643e7SJed Brown 
3316643e7SJed Brown /*
4316643e7SJed Brown   Code for timestepping with implicit Theta method
5316643e7SJed Brown 
6316643e7SJed Brown   Notes:
7316643e7SJed Brown   This method can be applied to DAE.
8316643e7SJed Brown 
9316643e7SJed Brown   This method is cast as a 1-stage implicit Runge-Kutta method.
10316643e7SJed Brown 
11316643e7SJed Brown   Theta | Theta
12316643e7SJed Brown   -------------
13316643e7SJed Brown         |  1
14316643e7SJed Brown 
15316643e7SJed Brown   To apply a diagonally implicit RK method to DAE, the stage formula
16316643e7SJed Brown 
17316643e7SJed Brown   X_i = x + h sum_j a_ij X'_j
18316643e7SJed Brown 
19316643e7SJed Brown   is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i)
20316643e7SJed Brown */
21316643e7SJed Brown #include "private/tsimpl.h"                /*I   "petscts.h"   I*/
22316643e7SJed Brown 
23316643e7SJed Brown typedef struct {
24316643e7SJed Brown   Vec X,Xdot;                   /* Storage for one stage */
25316643e7SJed Brown   Vec res;                      /* DAE residuals */
26ace68cafSJed Brown   PetscTruth extrapolate;
27316643e7SJed Brown   PetscReal Theta;
28316643e7SJed Brown   PetscReal shift;
29316643e7SJed Brown   PetscReal stage_time;
30316643e7SJed Brown } TS_Theta;
31316643e7SJed Brown 
32316643e7SJed Brown #undef __FUNCT__
33316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
34316643e7SJed Brown static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime)
35316643e7SJed Brown {
36316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
372b5a38e1SLisandro Dalcin   PetscInt       i,max_steps = ts->max_steps,its,lits;
382b5a38e1SLisandro Dalcin   PetscErrorCode ierr;
39316643e7SJed Brown 
40316643e7SJed Brown   PetscFunctionBegin;
41316643e7SJed Brown   *steps = -ts->steps;
422b5a38e1SLisandro Dalcin   *ptime = ts->ptime;
432b5a38e1SLisandro Dalcin 
442b5a38e1SLisandro Dalcin   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
45316643e7SJed Brown 
46316643e7SJed Brown   for (i=0; i<max_steps; i++) {
47316643e7SJed Brown     if (ts->ptime + ts->time_step > ts->max_time) break;
483f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
492b5a38e1SLisandro Dalcin 
50316643e7SJed Brown     th->stage_time = ts->ptime + th->Theta*ts->time_step;
51316643e7SJed Brown     th->shift = 1./(th->Theta*ts->time_step);
52316643e7SJed Brown 
53ace68cafSJed Brown     if (th->extrapolate) {
542b5a38e1SLisandro Dalcin       ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
55ace68cafSJed Brown     } else {
562b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
57ace68cafSJed Brown     }
58316643e7SJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr);
59316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
60316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
61316643e7SJed Brown     ts->nonlinear_its += its; ts->linear_its += lits;
622b5a38e1SLisandro Dalcin 
632b5a38e1SLisandro Dalcin     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
642b5a38e1SLisandro Dalcin     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
652b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
66316643e7SJed Brown     ts->steps++;
672b5a38e1SLisandro Dalcin 
683f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
692b5a38e1SLisandro Dalcin     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
70316643e7SJed Brown   }
71316643e7SJed Brown 
72316643e7SJed Brown   *steps += ts->steps;
73316643e7SJed Brown   *ptime  = ts->ptime;
74316643e7SJed Brown   PetscFunctionReturn(0);
75316643e7SJed Brown }
76316643e7SJed Brown 
77316643e7SJed Brown /*------------------------------------------------------------*/
78316643e7SJed Brown #undef __FUNCT__
79316643e7SJed Brown #define __FUNCT__ "TSDestroy_Theta"
80316643e7SJed Brown static PetscErrorCode TSDestroy_Theta(TS ts)
81316643e7SJed Brown {
82316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
83316643e7SJed Brown   PetscErrorCode  ierr;
84316643e7SJed Brown 
85316643e7SJed Brown   PetscFunctionBegin;
86a5cbb462SJed Brown   if (th->X)    {ierr = VecDestroy(th->X);CHKERRQ(ierr);}
87a5cbb462SJed Brown   if (th->Xdot) {ierr = VecDestroy(th->Xdot);CHKERRQ(ierr);}
88a5cbb462SJed Brown   if (th->res)  {ierr = VecDestroy(th->res);CHKERRQ(ierr);}
89316643e7SJed Brown   ierr = PetscFree(th);CHKERRQ(ierr);
90316643e7SJed Brown   PetscFunctionReturn(0);
91316643e7SJed Brown }
92316643e7SJed Brown 
93316643e7SJed Brown /*
94316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
952b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
96316643e7SJed Brown */
97316643e7SJed Brown #undef __FUNCT__
98*0f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
99*0f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
100316643e7SJed Brown {
101316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
102316643e7SJed Brown   PetscErrorCode ierr;
103316643e7SJed Brown 
104316643e7SJed Brown   PetscFunctionBegin;
1052b5a38e1SLisandro Dalcin   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr);
106316643e7SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr);
107316643e7SJed Brown   PetscFunctionReturn(0);
108316643e7SJed Brown }
109316643e7SJed Brown 
110316643e7SJed Brown #undef __FUNCT__
111*0f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
112*0f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
113316643e7SJed Brown {
114316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
115316643e7SJed Brown   PetscErrorCode ierr;
116316643e7SJed Brown 
117316643e7SJed Brown   PetscFunctionBegin;
118*0f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
119316643e7SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);CHKERRQ(ierr);
120316643e7SJed Brown   PetscFunctionReturn(0);
121316643e7SJed Brown }
122316643e7SJed Brown 
123316643e7SJed Brown 
124316643e7SJed Brown #undef __FUNCT__
125316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
126316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
127316643e7SJed Brown {
128316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
129316643e7SJed Brown   PetscErrorCode ierr;
130316643e7SJed Brown 
131316643e7SJed Brown   PetscFunctionBegin;
1322b5a38e1SLisandro Dalcin   if (ts->problem_type == TS_LINEAR) {
1332b5a38e1SLisandro Dalcin     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
1342b5a38e1SLisandro Dalcin   }
135316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
136316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
137316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->res);CHKERRQ(ierr);
138*0f5c6efeSJed Brown   ierr = SNESSetFunction(ts->snes,th->res,SNESTSFormFunction,ts);CHKERRQ(ierr);
1397c0b301bSJed Brown   /* This is nasty.  SNESSetFromOptions() is usually called in TSSetFromOptions().  With -snes_mf_operator, it will
1407c0b301bSJed Brown   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
1417c0b301bSJed Brown   context.  Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets
1427c0b301bSJed Brown   the Jacobian user context to snes->funP, it will actually be NULL.  This is not a problem because both snes->funP and
1437c0b301bSJed Brown   snes->jacP should be the TS. */
1447c0b301bSJed Brown   {
1457c0b301bSJed Brown     Mat A,B;
1467c0b301bSJed Brown     PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
1477c0b301bSJed Brown     void *ctx;
1487c0b301bSJed Brown     ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr);
149*0f5c6efeSJed Brown     ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr);
1507c0b301bSJed Brown   }
151316643e7SJed Brown   PetscFunctionReturn(0);
152316643e7SJed Brown }
153316643e7SJed Brown /*------------------------------------------------------------*/
154316643e7SJed Brown 
155316643e7SJed Brown #undef __FUNCT__
156316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
157316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
158316643e7SJed Brown {
159316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
160316643e7SJed Brown   PetscErrorCode ierr;
161316643e7SJed Brown 
162316643e7SJed Brown   PetscFunctionBegin;
163d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
164316643e7SJed Brown   {
165316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
166ace68cafSJed Brown     ierr = PetscOptionsTruth("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
167316643e7SJed Brown   }
168316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
169316643e7SJed Brown   PetscFunctionReturn(0);
170316643e7SJed Brown }
171316643e7SJed Brown 
172316643e7SJed Brown #undef __FUNCT__
173316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
174316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
175316643e7SJed Brown {
176316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
177316643e7SJed Brown   PetscTruth      iascii;
178316643e7SJed Brown   PetscErrorCode  ierr;
179316643e7SJed Brown 
180316643e7SJed Brown   PetscFunctionBegin;
181316643e7SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
182316643e7SJed Brown   if (iascii) {
183316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
184ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
185316643e7SJed Brown   } else {
186316643e7SJed Brown     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name);
187316643e7SJed Brown   }
188316643e7SJed Brown   PetscFunctionReturn(0);
189316643e7SJed Brown }
190316643e7SJed Brown 
191316643e7SJed Brown /* ------------------------------------------------------------ */
192316643e7SJed Brown /*MC
19396f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
194316643e7SJed Brown 
195316643e7SJed Brown   Level: beginner
196316643e7SJed Brown 
197316643e7SJed Brown .seealso:  TSCreate(), TS, TSSetType()
198316643e7SJed Brown 
199316643e7SJed Brown M*/
200316643e7SJed Brown EXTERN_C_BEGIN
201316643e7SJed Brown #undef __FUNCT__
202316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
203316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Theta(TS ts)
204316643e7SJed Brown {
205316643e7SJed Brown   TS_Theta       *th;
206316643e7SJed Brown   PetscErrorCode ierr;
207316643e7SJed Brown 
208316643e7SJed Brown   PetscFunctionBegin;
209316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
210316643e7SJed Brown   ts->ops->view           = TSView_Theta;
211316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
212316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
213316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
214*0f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
215*0f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
216316643e7SJed Brown 
2172b5a38e1SLisandro Dalcin   ts->problem_type = TS_NONLINEAR;
218316643e7SJed Brown   ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
219316643e7SJed Brown   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
220316643e7SJed Brown 
221316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
222316643e7SJed Brown   ts->data = (void*)th;
223316643e7SJed Brown 
224ace68cafSJed Brown   th->extrapolate = PETSC_TRUE;
225316643e7SJed Brown   th->Theta       = 0.5;
226316643e7SJed Brown 
227316643e7SJed Brown   PetscFunctionReturn(0);
228316643e7SJed Brown }
229316643e7SJed Brown EXTERN_C_END
230