xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision c6db04a5321582041def2b1e244c75985478b3ef)
1316643e7SJed Brown 
2316643e7SJed Brown /*
3316643e7SJed Brown   Code for timestepping with implicit Theta method
4316643e7SJed Brown 
5316643e7SJed Brown   Notes:
6316643e7SJed Brown   This method can be applied to DAE.
7316643e7SJed Brown 
8316643e7SJed Brown   This method is cast as a 1-stage implicit Runge-Kutta method.
9316643e7SJed Brown 
10316643e7SJed Brown   Theta | Theta
11316643e7SJed Brown   -------------
12316643e7SJed Brown         |  1
13316643e7SJed Brown 
14316643e7SJed Brown   To apply a diagonally implicit RK method to DAE, the stage formula
15316643e7SJed Brown 
16316643e7SJed Brown   X_i = x + h sum_j a_ij X'_j
17316643e7SJed Brown 
18316643e7SJed Brown   is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i)
19316643e7SJed Brown */
20*c6db04a5SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
21316643e7SJed Brown 
22316643e7SJed Brown typedef struct {
23316643e7SJed Brown   Vec X,Xdot;                   /* Storage for one stage */
24ace3abfcSBarry Smith   PetscBool  extrapolate;
25316643e7SJed Brown   PetscReal Theta;
26316643e7SJed Brown   PetscReal shift;
27316643e7SJed Brown   PetscReal stage_time;
28316643e7SJed Brown } TS_Theta;
29316643e7SJed Brown 
30316643e7SJed Brown #undef __FUNCT__
31316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
32316643e7SJed Brown static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime)
33316643e7SJed Brown {
34316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
352b5a38e1SLisandro Dalcin   PetscInt       i,max_steps = ts->max_steps,its,lits;
362b5a38e1SLisandro Dalcin   PetscErrorCode ierr;
37316643e7SJed Brown 
38316643e7SJed Brown   PetscFunctionBegin;
39316643e7SJed Brown   *steps = -ts->steps;
402b5a38e1SLisandro Dalcin   *ptime = ts->ptime;
412b5a38e1SLisandro Dalcin 
422b5a38e1SLisandro Dalcin   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
43316643e7SJed Brown 
44316643e7SJed Brown   for (i=0; i<max_steps; i++) {
45316643e7SJed Brown     if (ts->ptime + ts->time_step > ts->max_time) break;
463f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
472b5a38e1SLisandro Dalcin 
48316643e7SJed Brown     th->stage_time = ts->ptime + th->Theta*ts->time_step;
49316643e7SJed Brown     th->shift = 1./(th->Theta*ts->time_step);
50316643e7SJed Brown 
51ace68cafSJed Brown     if (th->extrapolate) {
522b5a38e1SLisandro Dalcin       ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
53ace68cafSJed Brown     } else {
542b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
55ace68cafSJed Brown     }
56316643e7SJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr);
57316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
58316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
59316643e7SJed Brown     ts->nonlinear_its += its; ts->linear_its += lits;
602b5a38e1SLisandro Dalcin 
612b5a38e1SLisandro Dalcin     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
622b5a38e1SLisandro Dalcin     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
632b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
64316643e7SJed Brown     ts->steps++;
652b5a38e1SLisandro Dalcin 
663f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
672b5a38e1SLisandro Dalcin     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
68316643e7SJed Brown   }
69316643e7SJed Brown 
70316643e7SJed Brown   *steps += ts->steps;
71316643e7SJed Brown   *ptime  = ts->ptime;
72316643e7SJed Brown   PetscFunctionReturn(0);
73316643e7SJed Brown }
74316643e7SJed Brown 
75316643e7SJed Brown /*------------------------------------------------------------*/
76316643e7SJed Brown #undef __FUNCT__
77316643e7SJed Brown #define __FUNCT__ "TSDestroy_Theta"
78316643e7SJed Brown static PetscErrorCode TSDestroy_Theta(TS ts)
79316643e7SJed Brown {
80316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
81316643e7SJed Brown   PetscErrorCode  ierr;
82316643e7SJed Brown 
83316643e7SJed Brown   PetscFunctionBegin;
84a5cbb462SJed Brown   if (th->X)    {ierr = VecDestroy(th->X);CHKERRQ(ierr);}
85a5cbb462SJed Brown   if (th->Xdot) {ierr = VecDestroy(th->Xdot);CHKERRQ(ierr);}
86316643e7SJed Brown   ierr = PetscFree(th);CHKERRQ(ierr);
87316643e7SJed Brown   PetscFunctionReturn(0);
88316643e7SJed Brown }
89316643e7SJed Brown 
90316643e7SJed Brown /*
91316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
922b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
93316643e7SJed Brown */
94316643e7SJed Brown #undef __FUNCT__
950f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
960f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
97316643e7SJed Brown {
98316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
99316643e7SJed Brown   PetscErrorCode ierr;
100316643e7SJed Brown 
101316643e7SJed Brown   PetscFunctionBegin;
1022b5a38e1SLisandro Dalcin   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr);
103316643e7SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr);
104316643e7SJed Brown   PetscFunctionReturn(0);
105316643e7SJed Brown }
106316643e7SJed Brown 
107316643e7SJed Brown #undef __FUNCT__
1080f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
1090f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
110316643e7SJed Brown {
111316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
112316643e7SJed Brown   PetscErrorCode ierr;
113316643e7SJed Brown 
114316643e7SJed Brown   PetscFunctionBegin;
1150f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
116316643e7SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);CHKERRQ(ierr);
117316643e7SJed Brown   PetscFunctionReturn(0);
118316643e7SJed Brown }
119316643e7SJed Brown 
120316643e7SJed Brown 
121316643e7SJed Brown #undef __FUNCT__
122316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
123316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
124316643e7SJed Brown {
125316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
126316643e7SJed Brown   PetscErrorCode ierr;
127a6e0575dSJed Brown   Vec            res;
128316643e7SJed Brown 
129316643e7SJed Brown   PetscFunctionBegin;
13017186662SBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
131316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
132316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
133a6e0575dSJed Brown   ierr = VecDuplicate(ts->vec_sol,&res);CHKERRQ(ierr);
134a6e0575dSJed Brown   ierr = SNESSetFunction(ts->snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
135a6e0575dSJed Brown   ierr = VecDestroy(res);CHKERRQ(ierr);
1367c0b301bSJed Brown   /* This is nasty.  SNESSetFromOptions() is usually called in TSSetFromOptions().  With -snes_mf_operator, it will
1377c0b301bSJed 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
1387c0b301bSJed Brown   context.  Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets
1397c0b301bSJed Brown   the Jacobian user context to snes->funP, it will actually be NULL.  This is not a problem because both snes->funP and
1407c0b301bSJed Brown   snes->jacP should be the TS. */
1417c0b301bSJed Brown   {
1427c0b301bSJed Brown     Mat A,B;
1437c0b301bSJed Brown     PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
1447c0b301bSJed Brown     void *ctx;
1457c0b301bSJed Brown     ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr);
1460f5c6efeSJed Brown     ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr);
1477c0b301bSJed Brown   }
148316643e7SJed Brown   PetscFunctionReturn(0);
149316643e7SJed Brown }
150316643e7SJed Brown /*------------------------------------------------------------*/
151316643e7SJed Brown 
152316643e7SJed Brown #undef __FUNCT__
153316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
154316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
155316643e7SJed Brown {
156316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
157316643e7SJed Brown   PetscErrorCode ierr;
158316643e7SJed Brown 
159316643e7SJed Brown   PetscFunctionBegin;
160d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
161316643e7SJed Brown   {
162316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
163acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
164316643e7SJed Brown   }
165316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
166316643e7SJed Brown   PetscFunctionReturn(0);
167316643e7SJed Brown }
168316643e7SJed Brown 
169316643e7SJed Brown #undef __FUNCT__
170316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
171316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
172316643e7SJed Brown {
173316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
174ace3abfcSBarry Smith   PetscBool       iascii;
175316643e7SJed Brown   PetscErrorCode  ierr;
176316643e7SJed Brown 
177316643e7SJed Brown   PetscFunctionBegin;
1782692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
179316643e7SJed Brown   if (iascii) {
180316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
181ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
182316643e7SJed Brown   } else {
183e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name);
184316643e7SJed Brown   }
185316643e7SJed Brown   PetscFunctionReturn(0);
186316643e7SJed Brown }
187316643e7SJed Brown 
1880de4c49aSJed Brown EXTERN_C_BEGIN
1890de4c49aSJed Brown #undef __FUNCT__
1900de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
1917087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
1920de4c49aSJed Brown {
1930de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1940de4c49aSJed Brown 
1950de4c49aSJed Brown   PetscFunctionBegin;
1960de4c49aSJed Brown   *theta = th->Theta;
1970de4c49aSJed Brown   PetscFunctionReturn(0);
1980de4c49aSJed Brown }
1990de4c49aSJed Brown 
2000de4c49aSJed Brown #undef __FUNCT__
2010de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
2027087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
2030de4c49aSJed Brown {
2040de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
2050de4c49aSJed Brown 
2060de4c49aSJed Brown   PetscFunctionBegin;
207e7be1afaSJed Brown   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
2080de4c49aSJed Brown   th->Theta = theta;
2090de4c49aSJed Brown   PetscFunctionReturn(0);
2100de4c49aSJed Brown }
2110de4c49aSJed Brown EXTERN_C_END
2120de4c49aSJed Brown 
213316643e7SJed Brown /* ------------------------------------------------------------ */
214316643e7SJed Brown /*MC
21596f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
216316643e7SJed Brown 
217316643e7SJed Brown   Level: beginner
218316643e7SJed Brown 
219316643e7SJed Brown .seealso:  TSCreate(), TS, TSSetType()
220316643e7SJed Brown 
221316643e7SJed Brown M*/
222316643e7SJed Brown EXTERN_C_BEGIN
223316643e7SJed Brown #undef __FUNCT__
224316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
2257087cfbeSBarry Smith PetscErrorCode  TSCreate_Theta(TS ts)
226316643e7SJed Brown {
227316643e7SJed Brown   TS_Theta       *th;
228316643e7SJed Brown   PetscErrorCode ierr;
229316643e7SJed Brown 
230316643e7SJed Brown   PetscFunctionBegin;
231316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
232316643e7SJed Brown   ts->ops->view           = TSView_Theta;
233316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
234316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
235316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
2360f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
2370f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
238316643e7SJed Brown 
2392b5a38e1SLisandro Dalcin   ts->problem_type = TS_NONLINEAR;
240316643e7SJed Brown   ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
241316643e7SJed Brown   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
242316643e7SJed Brown 
243316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
244316643e7SJed Brown   ts->data = (void*)th;
245316643e7SJed Brown 
2466f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
247316643e7SJed Brown   th->Theta       = 0.5;
248316643e7SJed Brown 
2490de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
2500de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
251316643e7SJed Brown   PetscFunctionReturn(0);
252316643e7SJed Brown }
253316643e7SJed Brown EXTERN_C_END
2540de4c49aSJed Brown 
2550de4c49aSJed Brown #undef __FUNCT__
2560de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
2570de4c49aSJed Brown /*@
2580de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
2590de4c49aSJed Brown 
2600de4c49aSJed Brown   Not Collective
2610de4c49aSJed Brown 
2620de4c49aSJed Brown   Input Parameter:
2630de4c49aSJed Brown .  ts - timestepping context
2640de4c49aSJed Brown 
2650de4c49aSJed Brown   Output Parameter:
2660de4c49aSJed Brown .  theta - stage abscissa
2670de4c49aSJed Brown 
2680de4c49aSJed Brown   Note:
2690de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
2700de4c49aSJed Brown 
2710de4c49aSJed Brown   Level: Advanced
2720de4c49aSJed Brown 
2730de4c49aSJed Brown .seealso: TSThetaSetTheta()
2740de4c49aSJed Brown @*/
2757087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
2760de4c49aSJed Brown {
2774ac538c5SBarry Smith   PetscErrorCode ierr;
2780de4c49aSJed Brown 
2790de4c49aSJed Brown   PetscFunctionBegin;
280afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2810de4c49aSJed Brown   PetscValidPointer(theta,2);
2824ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
2830de4c49aSJed Brown   PetscFunctionReturn(0);
2840de4c49aSJed Brown }
2850de4c49aSJed Brown 
2860de4c49aSJed Brown #undef __FUNCT__
2870de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
2880de4c49aSJed Brown /*@
2890de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
2900de4c49aSJed Brown 
2910de4c49aSJed Brown   Not Collective
2920de4c49aSJed Brown 
2930de4c49aSJed Brown   Input Parameter:
2940de4c49aSJed Brown +  ts - timestepping context
2950de4c49aSJed Brown -  theta - stage abscissa
2960de4c49aSJed Brown 
2970de4c49aSJed Brown   Options Database:
2980de4c49aSJed Brown .  -ts_theta_theta <theta>
2990de4c49aSJed Brown 
3000de4c49aSJed Brown   Level: Intermediate
3010de4c49aSJed Brown 
3020de4c49aSJed Brown .seealso: TSThetaGetTheta()
3030de4c49aSJed Brown @*/
3047087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
3050de4c49aSJed Brown {
3064ac538c5SBarry Smith   PetscErrorCode ierr;
3070de4c49aSJed Brown 
3080de4c49aSJed Brown   PetscFunctionBegin;
309afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3104ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
3110de4c49aSJed Brown   PetscFunctionReturn(0);
3120de4c49aSJed Brown }
313