xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision cdbf8f939cdfb1c797c4b7f2cbbd00be19935363)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4c6db04a5SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
5316643e7SJed Brown 
6316643e7SJed Brown typedef struct {
7316643e7SJed Brown   Vec       X,Xdot;                   /* Storage for one stage */
8eb284becSJed Brown   Vec       affine;                   /* Affine vector needed for residual at beginning of step */
9ace3abfcSBarry Smith   PetscBool extrapolate;
10eb284becSJed Brown   PetscBool endpoint;
11316643e7SJed Brown   PetscReal Theta;
12316643e7SJed Brown   PetscReal shift;
13316643e7SJed Brown   PetscReal stage_time;
14316643e7SJed Brown } TS_Theta;
15316643e7SJed Brown 
16316643e7SJed Brown #undef __FUNCT__
17316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
18193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
19316643e7SJed Brown {
20316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
21b70ae86eSJed Brown   PetscInt       its,lits;
22*cdbf8f93SLisandro Dalcin   PetscReal      next_time_step;
232b5a38e1SLisandro Dalcin   PetscErrorCode ierr;
24316643e7SJed Brown 
25316643e7SJed Brown   PetscFunctionBegin;
26*cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
27eb284becSJed Brown   th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
28316643e7SJed Brown   th->shift = 1./(th->Theta*ts->time_step);
29316643e7SJed Brown 
30eb284becSJed Brown   if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
31eb284becSJed Brown     ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
32eb284becSJed Brown     if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
33eb284becSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
34eb284becSJed Brown     ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
35eb284becSJed Brown   }
36ace68cafSJed Brown   if (th->extrapolate) {
372b5a38e1SLisandro Dalcin     ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
38ace68cafSJed Brown   } else {
392b5a38e1SLisandro Dalcin     ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
40ace68cafSJed Brown   }
41eb284becSJed Brown   ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
42316643e7SJed Brown   ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
43316643e7SJed Brown   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
44316643e7SJed Brown   ts->nonlinear_its += its; ts->linear_its += lits;
452b5a38e1SLisandro Dalcin 
46eb284becSJed Brown   if (th->endpoint) {
47eb284becSJed Brown     ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
48eb284becSJed Brown   } else {
492b5a38e1SLisandro Dalcin     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
502b5a38e1SLisandro Dalcin     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
51eb284becSJed Brown   }
522b5a38e1SLisandro Dalcin   ts->ptime += ts->time_step;
53*cdbf8f93SLisandro Dalcin   ts->time_step_prev = ts->time_step;
54*cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
55316643e7SJed Brown   ts->steps++;
56316643e7SJed Brown   PetscFunctionReturn(0);
57316643e7SJed Brown }
58316643e7SJed Brown 
59cd652676SJed Brown #undef __FUNCT__
60cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
61cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
62cd652676SJed Brown {
63cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
645a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
65cd652676SJed Brown   PetscErrorCode ierr;
66cd652676SJed Brown 
67cd652676SJed Brown   PetscFunctionBegin;
68a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
695a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
705a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
71cd652676SJed Brown   PetscFunctionReturn(0);
72cd652676SJed Brown }
73cd652676SJed Brown 
74316643e7SJed Brown /*------------------------------------------------------------*/
75316643e7SJed Brown #undef __FUNCT__
76277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
77277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
78316643e7SJed Brown {
79316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
80316643e7SJed Brown   PetscErrorCode  ierr;
81316643e7SJed Brown 
82316643e7SJed Brown   PetscFunctionBegin;
836bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
846bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
85eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
86277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
87277b19d0SLisandro Dalcin }
88277b19d0SLisandro Dalcin 
89277b19d0SLisandro Dalcin #undef __FUNCT__
90277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
91277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
92277b19d0SLisandro Dalcin {
93277b19d0SLisandro Dalcin   PetscErrorCode  ierr;
94277b19d0SLisandro Dalcin 
95277b19d0SLisandro Dalcin   PetscFunctionBegin;
96277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
97277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
98335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
99335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
100eb284becSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
101316643e7SJed Brown   PetscFunctionReturn(0);
102316643e7SJed Brown }
103316643e7SJed Brown 
104316643e7SJed Brown /*
105316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
1062b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
107316643e7SJed Brown */
108316643e7SJed Brown #undef __FUNCT__
1090f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
1100f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
111316643e7SJed Brown {
112316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
113316643e7SJed Brown   PetscErrorCode ierr;
114316643e7SJed Brown 
115316643e7SJed Brown   PetscFunctionBegin;
1165a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
1172b5a38e1SLisandro Dalcin   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr);
118214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
119316643e7SJed Brown   PetscFunctionReturn(0);
120316643e7SJed Brown }
121316643e7SJed Brown 
122316643e7SJed Brown #undef __FUNCT__
1230f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
1240f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
125316643e7SJed Brown {
126316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
127316643e7SJed Brown   PetscErrorCode ierr;
128316643e7SJed Brown 
129316643e7SJed Brown   PetscFunctionBegin;
1300f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
131214bc6a2SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
132316643e7SJed Brown   PetscFunctionReturn(0);
133316643e7SJed Brown }
134316643e7SJed Brown 
135316643e7SJed Brown 
136316643e7SJed Brown #undef __FUNCT__
137316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
138316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
139316643e7SJed Brown {
140316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
141316643e7SJed Brown   PetscErrorCode ierr;
142316643e7SJed Brown 
143316643e7SJed Brown   PetscFunctionBegin;
144316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
145316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
146316643e7SJed Brown   PetscFunctionReturn(0);
147316643e7SJed Brown }
148316643e7SJed Brown /*------------------------------------------------------------*/
149316643e7SJed Brown 
150316643e7SJed Brown #undef __FUNCT__
151316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
152316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
153316643e7SJed Brown {
154316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
155316643e7SJed Brown   PetscErrorCode ierr;
156316643e7SJed Brown 
157316643e7SJed Brown   PetscFunctionBegin;
158d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
159316643e7SJed Brown   {
160316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
161acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
162eb284becSJed Brown     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,PETSC_NULL);CHKERRQ(ierr);
163d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);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   }
183d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
184316643e7SJed Brown   PetscFunctionReturn(0);
185316643e7SJed Brown }
186316643e7SJed Brown 
1870de4c49aSJed Brown EXTERN_C_BEGIN
1880de4c49aSJed Brown #undef __FUNCT__
1890de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
1907087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
1910de4c49aSJed Brown {
1920de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1930de4c49aSJed Brown 
1940de4c49aSJed Brown   PetscFunctionBegin;
1950de4c49aSJed Brown   *theta = th->Theta;
1960de4c49aSJed Brown   PetscFunctionReturn(0);
1970de4c49aSJed Brown }
1980de4c49aSJed Brown 
1990de4c49aSJed Brown #undef __FUNCT__
2000de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
2017087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
2020de4c49aSJed Brown {
2030de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
2040de4c49aSJed Brown 
2050de4c49aSJed Brown   PetscFunctionBegin;
206e7be1afaSJed Brown   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
2070de4c49aSJed Brown   th->Theta = theta;
2080de4c49aSJed Brown   PetscFunctionReturn(0);
2090de4c49aSJed Brown }
210eb284becSJed Brown 
211eb284becSJed Brown #undef __FUNCT__
212eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint_Theta"
213eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
214eb284becSJed Brown {
215eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
216eb284becSJed Brown 
217eb284becSJed Brown   PetscFunctionBegin;
218eb284becSJed Brown   th->endpoint = flg;
219eb284becSJed Brown   PetscFunctionReturn(0);
220eb284becSJed Brown }
2210de4c49aSJed Brown EXTERN_C_END
2220de4c49aSJed Brown 
223316643e7SJed Brown /* ------------------------------------------------------------ */
224316643e7SJed Brown /*MC
22596f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
226316643e7SJed Brown 
227316643e7SJed Brown    Level: beginner
228316643e7SJed Brown 
229eb284becSJed Brown    Notes:
230eb284becSJed Brown    This method can be applied to DAE.
231eb284becSJed Brown 
232eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
233eb284becSJed Brown 
234eb284becSJed Brown .vb
235eb284becSJed Brown   Theta | Theta
236eb284becSJed Brown   -------------
237eb284becSJed Brown         |  1
238eb284becSJed Brown .ve
239eb284becSJed Brown 
240eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
241eb284becSJed Brown 
242eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
243eb284becSJed Brown 
244eb284becSJed Brown .vb
245eb284becSJed Brown   0 | 0         0
246eb284becSJed Brown   1 | 1-Theta   Theta
247eb284becSJed Brown   -------------------
248eb284becSJed Brown     | 1-Theta   Theta
249eb284becSJed Brown .ve
250eb284becSJed Brown 
251eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
252eb284becSJed Brown 
253eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
254eb284becSJed Brown 
255eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
256eb284becSJed Brown 
257eb284becSJed Brown    is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i)
258eb284becSJed Brown 
259eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
260316643e7SJed Brown 
261316643e7SJed Brown M*/
262316643e7SJed Brown EXTERN_C_BEGIN
263316643e7SJed Brown #undef __FUNCT__
264316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
2657087cfbeSBarry Smith PetscErrorCode  TSCreate_Theta(TS ts)
266316643e7SJed Brown {
267316643e7SJed Brown   TS_Theta       *th;
268316643e7SJed Brown   PetscErrorCode ierr;
269316643e7SJed Brown 
270316643e7SJed Brown   PetscFunctionBegin;
271277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
272316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
273316643e7SJed Brown   ts->ops->view           = TSView_Theta;
274316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
275316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
276cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
277316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
2780f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
2790f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
280316643e7SJed Brown 
281316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
282316643e7SJed Brown   ts->data = (void*)th;
283316643e7SJed Brown 
2846f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
285316643e7SJed Brown   th->Theta       = 0.5;
286316643e7SJed Brown 
2870de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
2880de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
289eb284becSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
290316643e7SJed Brown   PetscFunctionReturn(0);
291316643e7SJed Brown }
292316643e7SJed Brown EXTERN_C_END
2930de4c49aSJed Brown 
2940de4c49aSJed Brown #undef __FUNCT__
2950de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
2960de4c49aSJed Brown /*@
2970de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
2980de4c49aSJed Brown 
2990de4c49aSJed Brown   Not Collective
3000de4c49aSJed Brown 
3010de4c49aSJed Brown   Input Parameter:
3020de4c49aSJed Brown .  ts - timestepping context
3030de4c49aSJed Brown 
3040de4c49aSJed Brown   Output Parameter:
3050de4c49aSJed Brown .  theta - stage abscissa
3060de4c49aSJed Brown 
3070de4c49aSJed Brown   Note:
3080de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
3090de4c49aSJed Brown 
3100de4c49aSJed Brown   Level: Advanced
3110de4c49aSJed Brown 
3120de4c49aSJed Brown .seealso: TSThetaSetTheta()
3130de4c49aSJed Brown @*/
3147087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
3150de4c49aSJed Brown {
3164ac538c5SBarry Smith   PetscErrorCode ierr;
3170de4c49aSJed Brown 
3180de4c49aSJed Brown   PetscFunctionBegin;
319afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3200de4c49aSJed Brown   PetscValidPointer(theta,2);
3214ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
3220de4c49aSJed Brown   PetscFunctionReturn(0);
3230de4c49aSJed Brown }
3240de4c49aSJed Brown 
3250de4c49aSJed Brown #undef __FUNCT__
3260de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
3270de4c49aSJed Brown /*@
3280de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
3290de4c49aSJed Brown 
3300de4c49aSJed Brown   Not Collective
3310de4c49aSJed Brown 
3320de4c49aSJed Brown   Input Parameter:
3330de4c49aSJed Brown +  ts - timestepping context
3340de4c49aSJed Brown -  theta - stage abscissa
3350de4c49aSJed Brown 
3360de4c49aSJed Brown   Options Database:
3370de4c49aSJed Brown .  -ts_theta_theta <theta>
3380de4c49aSJed Brown 
3390de4c49aSJed Brown   Level: Intermediate
3400de4c49aSJed Brown 
3410de4c49aSJed Brown .seealso: TSThetaGetTheta()
3420de4c49aSJed Brown @*/
3437087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
3440de4c49aSJed Brown {
3454ac538c5SBarry Smith   PetscErrorCode ierr;
3460de4c49aSJed Brown 
3470de4c49aSJed Brown   PetscFunctionBegin;
348afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3494ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
3500de4c49aSJed Brown   PetscFunctionReturn(0);
3510de4c49aSJed Brown }
352f33bbcb6SJed Brown 
353eb284becSJed Brown #undef __FUNCT__
354eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
355eb284becSJed Brown /*@
356eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
357eb284becSJed Brown 
358eb284becSJed Brown   Not Collective
359eb284becSJed Brown 
360eb284becSJed Brown   Input Parameter:
361eb284becSJed Brown +  ts - timestepping context
362eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
363eb284becSJed Brown 
364eb284becSJed Brown   Options Database:
365eb284becSJed Brown .  -ts_theta_endpoint <flg>
366eb284becSJed Brown 
367eb284becSJed Brown   Level: Intermediate
368eb284becSJed Brown 
369eb284becSJed Brown .seealso: TSTHETA, TSCN
370eb284becSJed Brown @*/
371eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
372eb284becSJed Brown {
373eb284becSJed Brown   PetscErrorCode ierr;
374eb284becSJed Brown 
375eb284becSJed Brown   PetscFunctionBegin;
376eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
377eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
378eb284becSJed Brown   PetscFunctionReturn(0);
379eb284becSJed Brown }
380eb284becSJed Brown 
381f33bbcb6SJed Brown /*
382f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
383f33bbcb6SJed Brown  * The creation functions for these specializations are below.
384f33bbcb6SJed Brown  */
385f33bbcb6SJed Brown 
386f33bbcb6SJed Brown #undef __FUNCT__
387f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
388f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
389f33bbcb6SJed Brown {
390d52bd9f3SBarry Smith   PetscErrorCode ierr;
391d52bd9f3SBarry Smith 
392f33bbcb6SJed Brown   PetscFunctionBegin;
393d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
394f33bbcb6SJed Brown   PetscFunctionReturn(0);
395f33bbcb6SJed Brown }
396f33bbcb6SJed Brown 
397f33bbcb6SJed Brown /*MC
398f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
399f33bbcb6SJed Brown 
400f33bbcb6SJed Brown   Level: beginner
401f33bbcb6SJed Brown 
402f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
403f33bbcb6SJed Brown 
404f33bbcb6SJed Brown M*/
405f33bbcb6SJed Brown EXTERN_C_BEGIN
406f33bbcb6SJed Brown #undef __FUNCT__
407f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
408f33bbcb6SJed Brown PetscErrorCode  TSCreate_BEuler(TS ts)
409f33bbcb6SJed Brown {
410f33bbcb6SJed Brown   PetscErrorCode ierr;
411f33bbcb6SJed Brown 
412f33bbcb6SJed Brown   PetscFunctionBegin;
413f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
414f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
415f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
416f33bbcb6SJed Brown   PetscFunctionReturn(0);
417f33bbcb6SJed Brown }
418f33bbcb6SJed Brown EXTERN_C_END
419f33bbcb6SJed Brown 
420f33bbcb6SJed Brown #undef __FUNCT__
421f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
422f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
423f33bbcb6SJed Brown {
424d52bd9f3SBarry Smith   PetscErrorCode ierr;
425d52bd9f3SBarry Smith 
426f33bbcb6SJed Brown   PetscFunctionBegin;
427d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
428f33bbcb6SJed Brown   PetscFunctionReturn(0);
429f33bbcb6SJed Brown }
430f33bbcb6SJed Brown 
431f33bbcb6SJed Brown /*MC
432f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
433f33bbcb6SJed Brown 
434f33bbcb6SJed Brown   Level: beginner
435f33bbcb6SJed Brown 
436f33bbcb6SJed Brown   Notes:
4377cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
4387cf5af47SJed Brown 
4397cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
440f33bbcb6SJed Brown 
441f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
442f33bbcb6SJed Brown 
443f33bbcb6SJed Brown M*/
444f33bbcb6SJed Brown EXTERN_C_BEGIN
445f33bbcb6SJed Brown #undef __FUNCT__
446f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
447f33bbcb6SJed Brown PetscErrorCode  TSCreate_CN(TS ts)
448f33bbcb6SJed Brown {
449f33bbcb6SJed Brown   PetscErrorCode ierr;
450f33bbcb6SJed Brown 
451f33bbcb6SJed Brown   PetscFunctionBegin;
452f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
453f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
454eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
455f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
456f33bbcb6SJed Brown   PetscFunctionReturn(0);
457f33bbcb6SJed Brown }
458f33bbcb6SJed Brown EXTERN_C_END
459