xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 193ac0bc5128976c4ec2c1a67f3f9cb026b77f22)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown 
4316643e7SJed Brown   Notes:
5316643e7SJed Brown   This method can be applied to DAE.
6316643e7SJed Brown 
7316643e7SJed Brown   This method is cast as a 1-stage implicit Runge-Kutta method.
8316643e7SJed Brown 
9316643e7SJed Brown   Theta | Theta
10316643e7SJed Brown   -------------
11316643e7SJed Brown         |  1
12316643e7SJed Brown 
13316643e7SJed Brown   To apply a diagonally implicit RK method to DAE, the stage formula
14316643e7SJed Brown 
15316643e7SJed Brown   X_i = x + h sum_j a_ij X'_j
16316643e7SJed Brown 
17316643e7SJed Brown   is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i)
18316643e7SJed Brown */
19c6db04a5SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
20316643e7SJed Brown 
21316643e7SJed Brown typedef struct {
22316643e7SJed Brown   Vec       X,Xdot;                   /* Storage for one stage */
23ace3abfcSBarry Smith   PetscBool extrapolate;
24316643e7SJed Brown   PetscReal Theta;
25316643e7SJed Brown   PetscReal shift;
26316643e7SJed Brown   PetscReal stage_time;
27316643e7SJed Brown } TS_Theta;
28316643e7SJed Brown 
29316643e7SJed Brown #undef __FUNCT__
30316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
31*193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
32316643e7SJed Brown {
33316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
34186e87acSLisandro Dalcin   PetscInt       i,its,lits;
352b5a38e1SLisandro Dalcin   PetscErrorCode ierr;
36316643e7SJed Brown 
37316643e7SJed Brown   PetscFunctionBegin;
38*193ac0bcSJed Brown   ts->time_step = ts->next_time_step;
39316643e7SJed Brown   th->stage_time = ts->ptime + th->Theta*ts->time_step;
40316643e7SJed Brown   th->shift = 1./(th->Theta*ts->time_step);
41316643e7SJed Brown 
42ace68cafSJed Brown   if (th->extrapolate) {
432b5a38e1SLisandro Dalcin     ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
44ace68cafSJed Brown   } else {
452b5a38e1SLisandro Dalcin     ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
46ace68cafSJed Brown   }
47316643e7SJed Brown   ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr);
48316643e7SJed Brown   ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
49316643e7SJed Brown   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
50316643e7SJed Brown   ts->nonlinear_its += its; ts->linear_its += lits;
512b5a38e1SLisandro Dalcin 
522b5a38e1SLisandro Dalcin   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
532b5a38e1SLisandro Dalcin   ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
542b5a38e1SLisandro Dalcin   ts->ptime          += ts->time_step;
55*193ac0bcSJed Brown   ts->next_time_step  = ts->time_step;
56316643e7SJed Brown   ts->steps++;
57316643e7SJed Brown   PetscFunctionReturn(0);
58316643e7SJed Brown }
59316643e7SJed Brown 
60316643e7SJed Brown /*------------------------------------------------------------*/
61316643e7SJed Brown #undef __FUNCT__
62277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
63277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
64316643e7SJed Brown {
65316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
66316643e7SJed Brown   PetscErrorCode  ierr;
67316643e7SJed Brown 
68316643e7SJed Brown   PetscFunctionBegin;
696bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
706bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
71277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
72277b19d0SLisandro Dalcin }
73277b19d0SLisandro Dalcin 
74277b19d0SLisandro Dalcin #undef __FUNCT__
75277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
76277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
77277b19d0SLisandro Dalcin {
78277b19d0SLisandro Dalcin   PetscErrorCode  ierr;
79277b19d0SLisandro Dalcin 
80277b19d0SLisandro Dalcin   PetscFunctionBegin;
81277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
82277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
83316643e7SJed Brown   PetscFunctionReturn(0);
84316643e7SJed Brown }
85316643e7SJed Brown 
86316643e7SJed Brown /*
87316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
882b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
89316643e7SJed Brown */
90316643e7SJed Brown #undef __FUNCT__
910f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
920f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
93316643e7SJed Brown {
94316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
95316643e7SJed Brown   PetscErrorCode ierr;
96316643e7SJed Brown 
97316643e7SJed Brown   PetscFunctionBegin;
982b5a38e1SLisandro Dalcin   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr);
99214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
100316643e7SJed Brown   PetscFunctionReturn(0);
101316643e7SJed Brown }
102316643e7SJed Brown 
103316643e7SJed Brown #undef __FUNCT__
1040f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
1050f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
106316643e7SJed Brown {
107316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
108316643e7SJed Brown   PetscErrorCode ierr;
109316643e7SJed Brown 
110316643e7SJed Brown   PetscFunctionBegin;
1110f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
112214bc6a2SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
113316643e7SJed Brown   PetscFunctionReturn(0);
114316643e7SJed Brown }
115316643e7SJed Brown 
116316643e7SJed Brown 
117316643e7SJed Brown #undef __FUNCT__
118316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
119316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
120316643e7SJed Brown {
121316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
122316643e7SJed Brown   PetscErrorCode ierr;
123316643e7SJed Brown 
124316643e7SJed Brown   PetscFunctionBegin;
125316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
126316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
127316643e7SJed Brown   PetscFunctionReturn(0);
128316643e7SJed Brown }
129316643e7SJed Brown /*------------------------------------------------------------*/
130316643e7SJed Brown 
131316643e7SJed Brown #undef __FUNCT__
132316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
133316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
134316643e7SJed Brown {
135316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
136316643e7SJed Brown   PetscErrorCode ierr;
137316643e7SJed Brown 
138316643e7SJed Brown   PetscFunctionBegin;
139d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
140316643e7SJed Brown   {
141316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
142acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
143316643e7SJed Brown   }
144316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
145316643e7SJed Brown   PetscFunctionReturn(0);
146316643e7SJed Brown }
147316643e7SJed Brown 
148316643e7SJed Brown #undef __FUNCT__
149316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
150316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
151316643e7SJed Brown {
152316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
153ace3abfcSBarry Smith   PetscBool       iascii;
154316643e7SJed Brown   PetscErrorCode  ierr;
155316643e7SJed Brown 
156316643e7SJed Brown   PetscFunctionBegin;
1572692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
158316643e7SJed Brown   if (iascii) {
159316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
160ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
161316643e7SJed Brown   } else {
162e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name);
163316643e7SJed Brown   }
164316643e7SJed Brown   PetscFunctionReturn(0);
165316643e7SJed Brown }
166316643e7SJed Brown 
1670de4c49aSJed Brown EXTERN_C_BEGIN
1680de4c49aSJed Brown #undef __FUNCT__
1690de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
1707087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
1710de4c49aSJed Brown {
1720de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1730de4c49aSJed Brown 
1740de4c49aSJed Brown   PetscFunctionBegin;
1750de4c49aSJed Brown   *theta = th->Theta;
1760de4c49aSJed Brown   PetscFunctionReturn(0);
1770de4c49aSJed Brown }
1780de4c49aSJed Brown 
1790de4c49aSJed Brown #undef __FUNCT__
1800de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
1817087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
1820de4c49aSJed Brown {
1830de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1840de4c49aSJed Brown 
1850de4c49aSJed Brown   PetscFunctionBegin;
186e7be1afaSJed Brown   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
1870de4c49aSJed Brown   th->Theta = theta;
1880de4c49aSJed Brown   PetscFunctionReturn(0);
1890de4c49aSJed Brown }
1900de4c49aSJed Brown EXTERN_C_END
1910de4c49aSJed Brown 
192316643e7SJed Brown /* ------------------------------------------------------------ */
193316643e7SJed Brown /*MC
19496f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
195316643e7SJed Brown 
196316643e7SJed Brown   Level: beginner
197316643e7SJed Brown 
198316643e7SJed Brown .seealso:  TSCreate(), TS, TSSetType()
199316643e7SJed Brown 
200316643e7SJed Brown M*/
201316643e7SJed Brown EXTERN_C_BEGIN
202316643e7SJed Brown #undef __FUNCT__
203316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
2047087cfbeSBarry Smith PetscErrorCode  TSCreate_Theta(TS ts)
205316643e7SJed Brown {
206316643e7SJed Brown   TS_Theta       *th;
207316643e7SJed Brown   PetscErrorCode ierr;
208316643e7SJed Brown 
209316643e7SJed Brown   PetscFunctionBegin;
210277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
211316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
212316643e7SJed Brown   ts->ops->view           = TSView_Theta;
213316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
214316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
215316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
2160f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
2170f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
218316643e7SJed Brown 
219316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
220316643e7SJed Brown   ts->data = (void*)th;
221316643e7SJed Brown 
2226f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
223316643e7SJed Brown   th->Theta       = 0.5;
224316643e7SJed Brown 
2250de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
2260de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
227316643e7SJed Brown   PetscFunctionReturn(0);
228316643e7SJed Brown }
229316643e7SJed Brown EXTERN_C_END
2300de4c49aSJed Brown 
2310de4c49aSJed Brown #undef __FUNCT__
2320de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
2330de4c49aSJed Brown /*@
2340de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
2350de4c49aSJed Brown 
2360de4c49aSJed Brown   Not Collective
2370de4c49aSJed Brown 
2380de4c49aSJed Brown   Input Parameter:
2390de4c49aSJed Brown .  ts - timestepping context
2400de4c49aSJed Brown 
2410de4c49aSJed Brown   Output Parameter:
2420de4c49aSJed Brown .  theta - stage abscissa
2430de4c49aSJed Brown 
2440de4c49aSJed Brown   Note:
2450de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
2460de4c49aSJed Brown 
2470de4c49aSJed Brown   Level: Advanced
2480de4c49aSJed Brown 
2490de4c49aSJed Brown .seealso: TSThetaSetTheta()
2500de4c49aSJed Brown @*/
2517087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
2520de4c49aSJed Brown {
2534ac538c5SBarry Smith   PetscErrorCode ierr;
2540de4c49aSJed Brown 
2550de4c49aSJed Brown   PetscFunctionBegin;
256afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2570de4c49aSJed Brown   PetscValidPointer(theta,2);
2584ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
2590de4c49aSJed Brown   PetscFunctionReturn(0);
2600de4c49aSJed Brown }
2610de4c49aSJed Brown 
2620de4c49aSJed Brown #undef __FUNCT__
2630de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
2640de4c49aSJed Brown /*@
2650de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
2660de4c49aSJed Brown 
2670de4c49aSJed Brown   Not Collective
2680de4c49aSJed Brown 
2690de4c49aSJed Brown   Input Parameter:
2700de4c49aSJed Brown +  ts - timestepping context
2710de4c49aSJed Brown -  theta - stage abscissa
2720de4c49aSJed Brown 
2730de4c49aSJed Brown   Options Database:
2740de4c49aSJed Brown .  -ts_theta_theta <theta>
2750de4c49aSJed Brown 
2760de4c49aSJed Brown   Level: Intermediate
2770de4c49aSJed Brown 
2780de4c49aSJed Brown .seealso: TSThetaGetTheta()
2790de4c49aSJed Brown @*/
2807087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
2810de4c49aSJed Brown {
2824ac538c5SBarry Smith   PetscErrorCode ierr;
2830de4c49aSJed Brown 
2840de4c49aSJed Brown   PetscFunctionBegin;
285afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2864ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
2870de4c49aSJed Brown   PetscFunctionReturn(0);
2880de4c49aSJed Brown }
289