xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision acfcf0e5ba359b58e6a8a7af3f239cd7334278e8)
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 */
25ace3abfcSBarry Smith   PetscBool  extrapolate;
26316643e7SJed Brown   PetscReal Theta;
27316643e7SJed Brown   PetscReal shift;
28316643e7SJed Brown   PetscReal stage_time;
29316643e7SJed Brown } TS_Theta;
30316643e7SJed Brown 
31316643e7SJed Brown #undef __FUNCT__
32316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
33316643e7SJed Brown static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime)
34316643e7SJed Brown {
35316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
362b5a38e1SLisandro Dalcin   PetscInt       i,max_steps = ts->max_steps,its,lits;
372b5a38e1SLisandro Dalcin   PetscErrorCode ierr;
38316643e7SJed Brown 
39316643e7SJed Brown   PetscFunctionBegin;
40316643e7SJed Brown   *steps = -ts->steps;
412b5a38e1SLisandro Dalcin   *ptime = ts->ptime;
422b5a38e1SLisandro Dalcin 
432b5a38e1SLisandro Dalcin   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
44316643e7SJed Brown 
45316643e7SJed Brown   for (i=0; i<max_steps; i++) {
46316643e7SJed Brown     if (ts->ptime + ts->time_step > ts->max_time) break;
473f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
482b5a38e1SLisandro Dalcin 
49316643e7SJed Brown     th->stage_time = ts->ptime + th->Theta*ts->time_step;
50316643e7SJed Brown     th->shift = 1./(th->Theta*ts->time_step);
51316643e7SJed Brown 
52ace68cafSJed Brown     if (th->extrapolate) {
532b5a38e1SLisandro Dalcin       ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
54ace68cafSJed Brown     } else {
552b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
56ace68cafSJed Brown     }
57316643e7SJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr);
58316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
59316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
60316643e7SJed Brown     ts->nonlinear_its += its; ts->linear_its += lits;
612b5a38e1SLisandro Dalcin 
622b5a38e1SLisandro Dalcin     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
632b5a38e1SLisandro Dalcin     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
642b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
65316643e7SJed Brown     ts->steps++;
662b5a38e1SLisandro Dalcin 
673f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
682b5a38e1SLisandro Dalcin     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
69316643e7SJed Brown   }
70316643e7SJed Brown 
71316643e7SJed Brown   *steps += ts->steps;
72316643e7SJed Brown   *ptime  = ts->ptime;
73316643e7SJed Brown   PetscFunctionReturn(0);
74316643e7SJed Brown }
75316643e7SJed Brown 
76316643e7SJed Brown /*------------------------------------------------------------*/
77316643e7SJed Brown #undef __FUNCT__
78316643e7SJed Brown #define __FUNCT__ "TSDestroy_Theta"
79316643e7SJed Brown static PetscErrorCode TSDestroy_Theta(TS ts)
80316643e7SJed Brown {
81316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
82316643e7SJed Brown   PetscErrorCode  ierr;
83316643e7SJed Brown 
84316643e7SJed Brown   PetscFunctionBegin;
85a5cbb462SJed Brown   if (th->X)    {ierr = VecDestroy(th->X);CHKERRQ(ierr);}
86a5cbb462SJed Brown   if (th->Xdot) {ierr = VecDestroy(th->Xdot);CHKERRQ(ierr);}
87316643e7SJed Brown   ierr = PetscFree(th);CHKERRQ(ierr);
88316643e7SJed Brown   PetscFunctionReturn(0);
89316643e7SJed Brown }
90316643e7SJed Brown 
91316643e7SJed Brown /*
92316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
932b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
94316643e7SJed Brown */
95316643e7SJed Brown #undef __FUNCT__
960f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
970f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
98316643e7SJed Brown {
99316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
100316643e7SJed Brown   PetscErrorCode ierr;
101316643e7SJed Brown 
102316643e7SJed Brown   PetscFunctionBegin;
1032b5a38e1SLisandro Dalcin   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr);
104316643e7SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr);
105316643e7SJed Brown   PetscFunctionReturn(0);
106316643e7SJed Brown }
107316643e7SJed Brown 
108316643e7SJed Brown #undef __FUNCT__
1090f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
1100f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
111316643e7SJed Brown {
112316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
113316643e7SJed Brown   PetscErrorCode ierr;
114316643e7SJed Brown 
115316643e7SJed Brown   PetscFunctionBegin;
1160f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
117316643e7SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);CHKERRQ(ierr);
118316643e7SJed Brown   PetscFunctionReturn(0);
119316643e7SJed Brown }
120316643e7SJed Brown 
121316643e7SJed Brown 
122316643e7SJed Brown #undef __FUNCT__
123316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
124316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
125316643e7SJed Brown {
126316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
127316643e7SJed Brown   PetscErrorCode ierr;
128a6e0575dSJed Brown   Vec            res;
129316643e7SJed Brown 
130316643e7SJed Brown   PetscFunctionBegin;
13117186662SBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
132316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
133316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
134a6e0575dSJed Brown   ierr = VecDuplicate(ts->vec_sol,&res);CHKERRQ(ierr);
135a6e0575dSJed Brown   ierr = SNESSetFunction(ts->snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
136a6e0575dSJed Brown   ierr = VecDestroy(res);CHKERRQ(ierr);
1377c0b301bSJed Brown   /* This is nasty.  SNESSetFromOptions() is usually called in TSSetFromOptions().  With -snes_mf_operator, it will
1387c0b301bSJed 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
1397c0b301bSJed Brown   context.  Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets
1407c0b301bSJed Brown   the Jacobian user context to snes->funP, it will actually be NULL.  This is not a problem because both snes->funP and
1417c0b301bSJed Brown   snes->jacP should be the TS. */
1427c0b301bSJed Brown   {
1437c0b301bSJed Brown     Mat A,B;
1447c0b301bSJed Brown     PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
1457c0b301bSJed Brown     void *ctx;
1467c0b301bSJed Brown     ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr);
1470f5c6efeSJed Brown     ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr);
1487c0b301bSJed Brown   }
149316643e7SJed Brown   PetscFunctionReturn(0);
150316643e7SJed Brown }
151316643e7SJed Brown /*------------------------------------------------------------*/
152316643e7SJed Brown 
153316643e7SJed Brown #undef __FUNCT__
154316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
155316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
156316643e7SJed Brown {
157316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
158316643e7SJed Brown   PetscErrorCode ierr;
159316643e7SJed Brown 
160316643e7SJed Brown   PetscFunctionBegin;
161d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
162316643e7SJed Brown   {
163316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
164*acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
165316643e7SJed Brown   }
166316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
167316643e7SJed Brown   PetscFunctionReturn(0);
168316643e7SJed Brown }
169316643e7SJed Brown 
170316643e7SJed Brown #undef __FUNCT__
171316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
172316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
173316643e7SJed Brown {
174316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
175ace3abfcSBarry Smith   PetscBool       iascii;
176316643e7SJed Brown   PetscErrorCode  ierr;
177316643e7SJed Brown 
178316643e7SJed Brown   PetscFunctionBegin;
1792692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
180316643e7SJed Brown   if (iascii) {
181316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
182ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
183316643e7SJed Brown   } else {
184e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name);
185316643e7SJed Brown   }
186316643e7SJed Brown   PetscFunctionReturn(0);
187316643e7SJed Brown }
188316643e7SJed Brown 
1890de4c49aSJed Brown EXTERN_C_BEGIN
1900de4c49aSJed Brown #undef __FUNCT__
1910de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
1920de4c49aSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
1930de4c49aSJed Brown {
1940de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1950de4c49aSJed Brown 
1960de4c49aSJed Brown   PetscFunctionBegin;
1970de4c49aSJed Brown   *theta = th->Theta;
1980de4c49aSJed Brown   PetscFunctionReturn(0);
1990de4c49aSJed Brown }
2000de4c49aSJed Brown 
2010de4c49aSJed Brown #undef __FUNCT__
2020de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
2030de4c49aSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSThetaSetTheta_Theta(TS ts,PetscReal theta)
2040de4c49aSJed Brown {
2050de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
2060de4c49aSJed Brown 
2070de4c49aSJed Brown   PetscFunctionBegin;
208e7be1afaSJed Brown   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
2090de4c49aSJed Brown   th->Theta = theta;
2100de4c49aSJed Brown   PetscFunctionReturn(0);
2110de4c49aSJed Brown }
2120de4c49aSJed Brown EXTERN_C_END
2130de4c49aSJed Brown 
214316643e7SJed Brown /* ------------------------------------------------------------ */
215316643e7SJed Brown /*MC
21696f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
217316643e7SJed Brown 
218316643e7SJed Brown   Level: beginner
219316643e7SJed Brown 
220316643e7SJed Brown .seealso:  TSCreate(), TS, TSSetType()
221316643e7SJed Brown 
222316643e7SJed Brown M*/
223316643e7SJed Brown EXTERN_C_BEGIN
224316643e7SJed Brown #undef __FUNCT__
225316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
226316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Theta(TS ts)
227316643e7SJed Brown {
228316643e7SJed Brown   TS_Theta       *th;
229316643e7SJed Brown   PetscErrorCode ierr;
230316643e7SJed Brown 
231316643e7SJed Brown   PetscFunctionBegin;
232316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
233316643e7SJed Brown   ts->ops->view           = TSView_Theta;
234316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
235316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
236316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
2370f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
2380f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
239316643e7SJed Brown 
2402b5a38e1SLisandro Dalcin   ts->problem_type = TS_NONLINEAR;
241316643e7SJed Brown   ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
242316643e7SJed Brown   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
243316643e7SJed Brown 
244316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
245316643e7SJed Brown   ts->data = (void*)th;
246316643e7SJed Brown 
2476f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
248316643e7SJed Brown   th->Theta       = 0.5;
249316643e7SJed Brown 
2500de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
2510de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
252316643e7SJed Brown   PetscFunctionReturn(0);
253316643e7SJed Brown }
254316643e7SJed Brown EXTERN_C_END
2550de4c49aSJed Brown 
2560de4c49aSJed Brown #undef __FUNCT__
2570de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
2580de4c49aSJed Brown /*@
2590de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
2600de4c49aSJed Brown 
2610de4c49aSJed Brown   Not Collective
2620de4c49aSJed Brown 
2630de4c49aSJed Brown   Input Parameter:
2640de4c49aSJed Brown .  ts - timestepping context
2650de4c49aSJed Brown 
2660de4c49aSJed Brown   Output Parameter:
2670de4c49aSJed Brown .  theta - stage abscissa
2680de4c49aSJed Brown 
2690de4c49aSJed Brown   Note:
2700de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
2710de4c49aSJed Brown 
2720de4c49aSJed Brown   Level: Advanced
2730de4c49aSJed Brown 
2740de4c49aSJed Brown .seealso: TSThetaSetTheta()
2750de4c49aSJed Brown @*/
2760de4c49aSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSThetaGetTheta(TS ts,PetscReal *theta)
2770de4c49aSJed Brown {
2784ac538c5SBarry Smith   PetscErrorCode ierr;
2790de4c49aSJed Brown 
2800de4c49aSJed Brown   PetscFunctionBegin;
281afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2820de4c49aSJed Brown   PetscValidPointer(theta,2);
2834ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
2840de4c49aSJed Brown   PetscFunctionReturn(0);
2850de4c49aSJed Brown }
2860de4c49aSJed Brown 
2870de4c49aSJed Brown #undef __FUNCT__
2880de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
2890de4c49aSJed Brown /*@
2900de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
2910de4c49aSJed Brown 
2920de4c49aSJed Brown   Not Collective
2930de4c49aSJed Brown 
2940de4c49aSJed Brown   Input Parameter:
2950de4c49aSJed Brown +  ts - timestepping context
2960de4c49aSJed Brown -  theta - stage abscissa
2970de4c49aSJed Brown 
2980de4c49aSJed Brown   Options Database:
2990de4c49aSJed Brown .  -ts_theta_theta <theta>
3000de4c49aSJed Brown 
3010de4c49aSJed Brown   Level: Intermediate
3020de4c49aSJed Brown 
3030de4c49aSJed Brown .seealso: TSThetaGetTheta()
3040de4c49aSJed Brown @*/
3050de4c49aSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSThetaSetTheta(TS ts,PetscReal theta)
3060de4c49aSJed Brown {
3074ac538c5SBarry Smith   PetscErrorCode ierr;
3080de4c49aSJed Brown 
3090de4c49aSJed Brown   PetscFunctionBegin;
310afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3114ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
3120de4c49aSJed Brown   PetscFunctionReturn(0);
3130de4c49aSJed Brown }
314