xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 7c0b301b739d24f219e3424569c5fc0502872d69)
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 Xold;
25316643e7SJed Brown   Vec X,Xdot;                   /* Storage for one stage */
26316643e7SJed Brown   Vec res;                      /* DAE residuals */
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   Vec            sol = ts->vec_sol;
37316643e7SJed Brown   PetscErrorCode ierr;
38316643e7SJed Brown   PetscInt       i,max_steps = ts->max_steps,its,lits;
39316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
40316643e7SJed Brown 
41316643e7SJed Brown   PetscFunctionBegin;
42316643e7SJed Brown   *steps = -ts->steps;
43316643e7SJed Brown   ierr = TSMonitor(ts,ts->steps,ts->ptime,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;
47316643e7SJed Brown     th->stage_time = ts->ptime + th->Theta*ts->time_step;
48316643e7SJed Brown     th->shift = 1./(th->Theta*ts->time_step);
49316643e7SJed Brown     ts->ptime += ts->time_step;
50316643e7SJed Brown 
51316643e7SJed Brown     ierr = VecCopy(sol,th->Xold);CHKERRQ(ierr); /* Used within function evalutaion */
52316643e7SJed Brown     ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,sol);CHKERRQ(ierr);
53316643e7SJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr);
54316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
55316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
56316643e7SJed Brown     ts->nonlinear_its += its; ts->linear_its += lits;
57316643e7SJed Brown     ierr = VecAXPY(sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
58316643e7SJed Brown     ts->steps++;
59316643e7SJed Brown     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
60316643e7SJed Brown   }
61316643e7SJed Brown 
62316643e7SJed Brown   *steps += ts->steps;
63316643e7SJed Brown   *ptime  = ts->ptime;
64316643e7SJed Brown   PetscFunctionReturn(0);
65316643e7SJed Brown }
66316643e7SJed Brown 
67316643e7SJed Brown /*------------------------------------------------------------*/
68316643e7SJed Brown #undef __FUNCT__
69316643e7SJed Brown #define __FUNCT__ "TSDestroy_Theta"
70316643e7SJed Brown static PetscErrorCode TSDestroy_Theta(TS ts)
71316643e7SJed Brown {
72316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
73316643e7SJed Brown   PetscErrorCode  ierr;
74316643e7SJed Brown 
75316643e7SJed Brown   PetscFunctionBegin;
76316643e7SJed Brown   ierr = VecDestroy(th->Xold);CHKERRQ(ierr);
77316643e7SJed Brown   ierr = VecDestroy(th->X);CHKERRQ(ierr);
78316643e7SJed Brown   ierr = VecDestroy(th->Xdot);CHKERRQ(ierr);
79316643e7SJed Brown   ierr = VecDestroy(th->res);CHKERRQ(ierr);
80316643e7SJed Brown   ierr = PetscFree(th);CHKERRQ(ierr);
81316643e7SJed Brown   PetscFunctionReturn(0);
82316643e7SJed Brown }
83316643e7SJed Brown 
84316643e7SJed Brown /*
85316643e7SJed Brown     This defines the nonlinear equation that is to be solved with SNES
86316643e7SJed Brown     G(U) = F[t0+T*dt, U, (U-U0)*shift] = 0
87316643e7SJed Brown */
88316643e7SJed Brown #undef __FUNCT__
89316643e7SJed Brown #define __FUNCT__ "TSThetaFunction"
90316643e7SJed Brown static PetscErrorCode TSThetaFunction(SNES snes,Vec x,Vec y,void *ctx)
91316643e7SJed Brown {
92316643e7SJed Brown   TS        ts = (TS)ctx;
93316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
94316643e7SJed Brown   PetscErrorCode ierr;
95316643e7SJed Brown 
96316643e7SJed Brown   PetscFunctionBegin;
97316643e7SJed Brown   ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->Xold,x);CHKERRQ(ierr);
98316643e7SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr);
99316643e7SJed Brown   PetscFunctionReturn(0);
100316643e7SJed Brown }
101316643e7SJed Brown 
102316643e7SJed Brown #undef __FUNCT__
103316643e7SJed Brown #define __FUNCT__ "TSThetaJacobian"
104316643e7SJed Brown static PetscErrorCode TSThetaJacobian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,void *ctx)
105316643e7SJed Brown {
106316643e7SJed Brown   TS        ts = (TS)ctx;
107316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
108316643e7SJed Brown   PetscErrorCode ierr;
109316643e7SJed Brown 
110316643e7SJed Brown   PetscFunctionBegin;
111316643e7SJed Brown   /* th->Xdot will have already been computed in TSThetaFunction */
112316643e7SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);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->Xold);CHKERRQ(ierr);
126316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
127316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
128316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->res);CHKERRQ(ierr);
129316643e7SJed Brown   ierr = SNESSetFunction(ts->snes,th->res,TSThetaFunction,ts);CHKERRQ(ierr);
130*7c0b301bSJed Brown   /* This is nasty.  SNESSetFromOptions() is usually called in TSSetFromOptions().  With -snes_mf_operator, it will
131*7c0b301bSJed 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
132*7c0b301bSJed Brown   context.  Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets
133*7c0b301bSJed Brown   the Jacobian user context to snes->funP, it will actually be NULL.  This is not a problem because both snes->funP and
134*7c0b301bSJed Brown   snes->jacP should be the TS. */
135*7c0b301bSJed Brown   {
136*7c0b301bSJed Brown     Mat A,B;
137*7c0b301bSJed Brown     PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
138*7c0b301bSJed Brown     void *ctx;
139*7c0b301bSJed Brown     ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr);
140*7c0b301bSJed Brown     ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&TSThetaJacobian,ctx?ctx:ts);CHKERRQ(ierr);
141*7c0b301bSJed Brown   }
142316643e7SJed Brown   PetscFunctionReturn(0);
143316643e7SJed Brown }
144316643e7SJed Brown /*------------------------------------------------------------*/
145316643e7SJed Brown 
146316643e7SJed Brown #undef __FUNCT__
147316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
148316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
149316643e7SJed Brown {
150316643e7SJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
151316643e7SJed Brown   PetscErrorCode ierr;
152316643e7SJed Brown 
153316643e7SJed Brown   PetscFunctionBegin;
154316643e7SJed Brown   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
155316643e7SJed Brown   {
156316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
157316643e7SJed Brown   }
158316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
159316643e7SJed Brown   PetscFunctionReturn(0);
160316643e7SJed Brown }
161316643e7SJed Brown 
162316643e7SJed Brown #undef __FUNCT__
163316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
164316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
165316643e7SJed Brown {
166316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
167316643e7SJed Brown   PetscTruth      iascii;
168316643e7SJed Brown   PetscErrorCode  ierr;
169316643e7SJed Brown 
170316643e7SJed Brown   PetscFunctionBegin;
171316643e7SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
172316643e7SJed Brown   if (iascii) {
173316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
174316643e7SJed Brown   } else {
175316643e7SJed Brown     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name);
176316643e7SJed Brown   }
177316643e7SJed Brown   PetscFunctionReturn(0);
178316643e7SJed Brown }
179316643e7SJed Brown 
180316643e7SJed Brown /* ------------------------------------------------------------ */
181316643e7SJed Brown /*MC
182*7c0b301bSJed Brown       TS_THETA - DAE solver using the implicit Theta method
183316643e7SJed Brown 
184316643e7SJed Brown   Level: beginner
185316643e7SJed Brown 
186316643e7SJed Brown .seealso:  TSCreate(), TS, TSSetType()
187316643e7SJed Brown 
188316643e7SJed Brown M*/
189316643e7SJed Brown EXTERN_C_BEGIN
190316643e7SJed Brown #undef __FUNCT__
191316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
192316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Theta(TS ts)
193316643e7SJed Brown {
194316643e7SJed Brown   TS_Theta       *th;
195316643e7SJed Brown   PetscErrorCode ierr;
196316643e7SJed Brown 
197316643e7SJed Brown   PetscFunctionBegin;
198316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
199316643e7SJed Brown   ts->ops->view           = TSView_Theta;
200316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
201316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
202316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
203316643e7SJed Brown 
204316643e7SJed Brown   ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
205316643e7SJed Brown   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
206316643e7SJed Brown 
207316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
208316643e7SJed Brown   ts->data = (void*)th;
209316643e7SJed Brown 
210316643e7SJed Brown   th->Theta = 0.5;
211316643e7SJed Brown 
212316643e7SJed Brown   PetscFunctionReturn(0);
213316643e7SJed Brown }
214316643e7SJed Brown EXTERN_C_END
215