xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 2e7b7f96234fe412e30f42a61403d5ed8380abcc)
1 /*
2   Code for timestepping with implicit Theta method
3 */
4 #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
5 #include <petscsnes.h>
6 #include <petscdm.h>
7 #include <petscmat.h>
8 
9 typedef struct {
10   /* context for time stepping */
11   PetscReal    stage_time;
12   Vec          X0,X,Xdot;                /* Storage for stages and time derivative */
13   Vec          affine;                   /* Affine vector needed for residual at beginning of step in endpoint formulation */
14   PetscReal    Theta;
15   PetscReal    ptime;
16   PetscReal    time_step;
17   PetscInt     order;
18   PetscBool    endpoint;
19   PetscBool    extrapolate;
20   TSStepStatus status;
21   Vec          VecCostIntegral0;         /* Backup for roll-backs due to events */
22 
23   /* context for sensitivity analysis */
24   PetscInt     num_tlm;                  /* Total number of tangent linear equations */
25   Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
26   Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
27   Vec          *VecsSensiTemp;           /* Vector to be multiplied with Jacobian transpose */
28   Mat          MatDeltaFwdSensip;        /* Increment of the forward sensitivity at stage */
29   Vec          VecDeltaFwdSensipCol;     /* Working vector for holding one column of the sensitivity matrix */
30   Mat          MatFwdSensip0;            /* backup for roll-backs due to events */
31   Vec          VecIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
32   Vec          *VecsIntegralSensip0;     /* backup for roll-backs due to events */
33   Vec          *VecsDeltaLam2;           /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
34   Vec          *VecsDeltaMu2;            /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
35   Vec          *VecsSensi2Temp;          /* Working vectors that holds the residual for the second-order adjoint */
36   Vec          *VecsAffine;              /* Working vectors to store residuals */
37   /* context for error estimation */
38   Vec          vec_sol_prev;
39   Vec          vec_lte_work;
40 } TS_Theta;
41 
42 static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
43 {
44   TS_Theta       *th = (TS_Theta*)ts->data;
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   if (X0) {
49     if (dm && dm != ts->dm) {
50       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
51     } else *X0 = ts->vec_sol;
52   }
53   if (Xdot) {
54     if (dm && dm != ts->dm) {
55       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
56     } else *Xdot = th->Xdot;
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
62 {
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   if (X0) {
67     if (dm && dm != ts->dm) {
68       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
69     }
70   }
71   if (Xdot) {
72     if (dm && dm != ts->dm) {
73       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
74     }
75   }
76   PetscFunctionReturn(0);
77 }
78 
79 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
80 {
81   PetscFunctionBegin;
82   PetscFunctionReturn(0);
83 }
84 
85 static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
86 {
87   TS             ts = (TS)ctx;
88   PetscErrorCode ierr;
89   Vec            X0,Xdot,X0_c,Xdot_c;
90 
91   PetscFunctionBegin;
92   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
93   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
94   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
95   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
96   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
97   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
98   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
99   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
100   PetscFunctionReturn(0);
101 }
102 
103 static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
104 {
105   PetscFunctionBegin;
106   PetscFunctionReturn(0);
107 }
108 
109 static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
110 {
111   TS             ts = (TS)ctx;
112   PetscErrorCode ierr;
113   Vec            X0,Xdot,X0_sub,Xdot_sub;
114 
115   PetscFunctionBegin;
116   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
117   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
118 
119   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121 
122   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124 
125   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
126   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
127   PetscFunctionReturn(0);
128 }
129 
130 static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
131 {
132   TS_Theta       *th = (TS_Theta*)ts->data;
133   PetscErrorCode ierr;
134 
135   PetscFunctionBegin;
136   if (th->endpoint) {
137     /* Evolve ts->vec_costintegral to compute integrals */
138     if (th->Theta!=1.0) {
139       ierr = TSComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
140       ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
141     }
142     ierr = TSComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
143     ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
144   } else {
145     ierr = TSComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
146     ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
147   }
148   PetscFunctionReturn(0);
149 }
150 
151 static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
152 {
153   TS_Theta       *th = (TS_Theta*)ts->data;
154   PetscErrorCode ierr;
155 
156   PetscFunctionBegin;
157   /* backup cost integral */
158   ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr);
159   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
160   PetscFunctionReturn(0);
161 }
162 
163 static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
164 {
165   PetscErrorCode ierr;
166 
167   PetscFunctionBegin;
168   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
169   PetscFunctionReturn(0);
170 }
171 
172 static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x)
173 {
174   PetscInt       nits,lits;
175   PetscErrorCode ierr;
176 
177   PetscFunctionBegin;
178   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
179   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
180   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
181   ts->snes_its += nits; ts->ksp_its += lits;
182   PetscFunctionReturn(0);
183 }
184 
185 static PetscErrorCode TSStep_Theta(TS ts)
186 {
187   TS_Theta       *th = (TS_Theta*)ts->data;
188   PetscInt       rejections = 0;
189   PetscBool      stageok,accept = PETSC_TRUE;
190   PetscReal      next_time_step = ts->time_step;
191   PetscErrorCode ierr;
192 
193   PetscFunctionBegin;
194   if (!ts->steprollback) {
195     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
196     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
197   }
198 
199   th->status = TS_STEP_INCOMPLETE;
200   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
201 
202     PetscReal shift = 1/(th->Theta*ts->time_step);
203     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
204 
205     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
206     if (th->extrapolate && !ts->steprestart) {
207       ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr);
208     }
209     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
210       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
211       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
212       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
213       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
214     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
215       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
216     }
217     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
218     ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
219     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
220     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
221     if (!stageok) goto reject_step;
222 
223     th->status = TS_STEP_PENDING;
224     if (th->endpoint) {
225       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
226     } else {
227       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr);
228       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
229     }
230     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
231     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
232     if (!accept) {
233       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
234       ts->time_step = next_time_step;
235       goto reject_step;
236     }
237 
238     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
239       th->ptime     = ts->ptime;
240       th->time_step = ts->time_step;
241     }
242 
243     ts->ptime += ts->time_step;
244     ts->time_step = next_time_step;
245     break;
246 
247   reject_step:
248     ts->reject++; accept = PETSC_FALSE;
249     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
250       ts->reason = TS_DIVERGED_STEP_REJECTED;
251       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
252     }
253   }
254   PetscFunctionReturn(0);
255 }
256 
257 static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
258 {
259   TS_Theta       *th = (TS_Theta*)ts->data;
260   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
261   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
262   PetscInt       nadj;
263   Mat            J,Jp;
264   KSP            ksp;
265   PetscReal      shift;
266   PetscScalar    *xarr;
267   PetscErrorCode ierr;
268 
269   PetscFunctionBegin;
270   th->status = TS_STEP_INCOMPLETE;
271   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
272   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
273 
274   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
275   th->stage_time = ts->ptime; /* time_step is negative*/
276   th->ptime      = ts->ptime + ts->time_step;
277   th->time_step  = -ts->time_step;
278 
279   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
280   ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr);
281   for (nadj=0; nadj<ts->numcost; nadj++) {
282     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
283     ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */
284     if (ts->vecs_drdu) {
285       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr);
286     }
287   }
288 
289   /* Build LHS for first-order adjoint */
290   shift = 1./th->time_step;
291   ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
292   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
293 
294   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
295   for (nadj=0; nadj<ts->numcost; nadj++) {
296     KSPConvergedReason kspreason;
297     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
298     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
299     if (kspreason < 0) {
300       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
301       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
302     }
303   }
304 
305   if (ts->vecs_sensi2) { /* U_{n+1} */
306     /* Get w1 at t_{n+1} from TLM matrix */
307     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
308     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
309     /* lambda_s^T F_UU w_1 */
310     ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
311     /* lambda_s^T F_UP w_2 */
312     ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
313     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
314       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
315       ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr);
316       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
317       if (ts->vecs_fup) {
318         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
319       }
320     }
321     /* Solve stage equation LHS X = RHS for second-order adjoint */
322     for (nadj=0; nadj<ts->numcost; nadj++) {
323       KSPConvergedReason kspreason;
324       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
325       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
326       if (kspreason < 0) {
327         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
328         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
329       }
330     }
331   }
332 
333   /* Update sensitivities, and evaluate integrals if there is any */
334   shift = 0.0;
335   ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_U */
336   ierr  = MatScale(J,-1.);CHKERRQ(ierr);
337   for (nadj=0; nadj<ts->numcost; nadj++) {
338     /* Add f_U \lambda_s to the original RHS */
339     ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
340     ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr);
341     ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
342     if (ts->vecs_sensi2) {
343       ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
344       ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr);
345       ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
346     }
347   }
348   if (ts->vecs_sensip) {
349     ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */
350     ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
351     if (ts->vecs_sensi2p) {
352       /* lambda_s^T F_PU w_1 */
353       ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
354       /* lambda_s^T F_PP w_2 */
355       ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
356     }
357     for (nadj=0; nadj<ts->numcost; nadj++) {
358       ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
359       ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
360       if (ts->vecs_drdp) {
361         ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
362       }
363       if (ts->vecs_sensi2p) {
364         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
365         ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
366         if (ts->vecs_fpu) {
367           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
368         }
369         if (ts->vecs_fpp) {
370           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
371         }
372       }
373     }
374   }
375 
376   if (ts->vecs_sensi2) {
377     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
378     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
379   }
380   th->status = TS_STEP_COMPLETE;
381   PetscFunctionReturn(0);
382 }
383 
384 static PetscErrorCode TSAdjointStep_Theta(TS ts)
385 {
386   TS_Theta       *th = (TS_Theta*)ts->data;
387   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
388   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
389   PetscInt       nadj;
390   Mat            J,Jp;
391   KSP            ksp;
392   PetscReal      shift;
393   PetscScalar    *xarr;
394   PetscErrorCode ierr;
395 
396   PetscFunctionBegin;
397   if (th->Theta == 1.) {
398     ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr);
399     PetscFunctionReturn(0);
400   }
401   th->status = TS_STEP_INCOMPLETE;
402   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
403   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
404 
405   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
406   th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/
407   th->ptime      = ts->ptime + ts->time_step;
408   th->time_step  = -ts->time_step;
409 
410   /* Build RHS for first-order adjoint */
411   /* Cost function has an integral term */
412   if (th->endpoint) {
413     ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr);
414   } else {
415     ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr);
416   }
417 
418   for (nadj=0; nadj<ts->numcost; nadj++) {
419     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
420     ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr);
421     if (ts->vec_costintegral) {
422       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr);
423     }
424   }
425 
426   /* Build LHS for first-order adjoint */
427   shift = 1./(th->Theta*th->time_step);
428   if (th->endpoint) {
429     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
430   } else {
431     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
432   }
433   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
434 
435   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
436   for (nadj=0; nadj<ts->numcost; nadj++) {
437     KSPConvergedReason kspreason;
438     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
439     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
440     if (kspreason < 0) {
441       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
442       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
443     }
444   }
445 
446   /* Second-order adjoint */
447   if (ts->vecs_sensi2) { /* U_{n+1} */
448     if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta");
449     /* Get w1 at t_{n+1} from TLM matrix */
450     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
451     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
452     /* lambda_s^T F_UU w_1 */
453     ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
454     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
455     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
456     /* lambda_s^T F_UP w_2 */
457     ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
458     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
459       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
460       ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr);
461       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
462       if (ts->vecs_fup) {
463         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
464       }
465     }
466     /* Solve stage equation LHS X = RHS for second-order adjoint */
467     for (nadj=0; nadj<ts->numcost; nadj++) {
468       KSPConvergedReason kspreason;
469       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
470       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
471       if (kspreason < 0) {
472         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
473         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
474       }
475     }
476   }
477 
478   /* Update sensitivities, and evaluate integrals if there is any */
479   if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
480     shift = 1./((th->Theta-1.)*th->time_step);
481     ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
482     /* R_U at t_n */
483     ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr);
484     for (nadj=0; nadj<ts->numcost; nadj++) {
485       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
486       if (ts->vecs_drdu) {
487         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr);
488       }
489       ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
490     }
491 
492     /* Second-order adjoint */
493     if (ts->vecs_sensi2) { /* U_n */
494       /* Get w1 at t_n from TLM matrix */
495       ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
496       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
497       /* lambda_s^T F_UU w_1 */
498       ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
499       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
500       ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
501       /* lambda_s^T F_UU w_2 */
502       ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
503       for (nadj=0; nadj<ts->numcost; nadj++) {
504         /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2  */
505         ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
506         ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
507         if (ts->vecs_fup) {
508           ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
509         }
510         ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr);
511       }
512     }
513 
514     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
515       /* U_{n+1} */
516       shift = -1./(th->Theta*th->time_step);
517       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
518       ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
519       for (nadj=0; nadj<ts->numcost; nadj++) {
520         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
521         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
522       }
523       if (ts->vecs_sensi2p) { /* second-order */
524         /* Get w1 at t_{n+1} from TLM matrix */
525         ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
526         ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
527         /* lambda_s^T F_PU w_1 */
528         ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
529         ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
530         ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
531 
532         /* lambda_s^T F_PP w_2 */
533         ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
534         for (nadj=0; nadj<ts->numcost; nadj++) {
535           /* Mu2 <- Mu2 + h theta F_P^T Lambda_s + h theta (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2)  */
536           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
537           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
538           if (ts->vecs_fpu) {
539             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
540           }
541           if (ts->vecs_fpp) {
542             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
543           }
544         }
545       }
546 
547       /* U_s */
548       shift = 1./((th->Theta-1.0)*th->time_step);
549       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
550       ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
551       for (nadj=0; nadj<ts->numcost; nadj++) {
552         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
553         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
554         if (ts->vecs_sensi2p) { /* second-order */
555           /* Get w1 at t_n from TLM matrix */
556           ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
557           ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
558           /* lambda_s^T F_PU w_1 */
559           ierr = TSComputeIHessianProductFunction3(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
560           ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
561           ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
562           /* lambda_s^T F_PP w_2 */
563           ierr = TSComputeIHessianProductFunction4(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
564           for (nadj=0; nadj<ts->numcost; nadj++) {
565             /* Mu2 <- Mu2 + h(1-theta) F_P^T Lambda_s + h(1-theta) (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */
566             ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
567             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr);
568             if (ts->vecs_fpu) {
569               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr);
570             }
571             if (ts->vecs_fpp) {
572               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr);
573             }
574           }
575         }
576       }
577     }
578   } else { /* one-stage case */
579     shift = 0.0;
580     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
581     ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr);
582     for (nadj=0; nadj<ts->numcost; nadj++) {
583       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
584       ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
585       if (ts->vecs_drdu) {
586         ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr);
587       }
588     }
589     if (ts->vecs_sensip) {
590       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
591       ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
592       for (nadj=0; nadj<ts->numcost; nadj++) {
593         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
594         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
595         if (ts->vecs_drdp) {
596           ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
597         }
598       }
599     }
600   }
601 
602   th->status = TS_STEP_COMPLETE;
603   PetscFunctionReturn(0);
604 }
605 
606 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
607 {
608   TS_Theta       *th = (TS_Theta*)ts->data;
609   PetscReal      dt  = t - ts->ptime;
610   PetscErrorCode ierr;
611 
612   PetscFunctionBegin;
613   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
614   if (th->endpoint) dt *= th->Theta;
615   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
616   PetscFunctionReturn(0);
617 }
618 
619 static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
620 {
621   TS_Theta       *th = (TS_Theta*)ts->data;
622   Vec            X = ts->vec_sol;      /* X = solution */
623   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
624   PetscReal      wltea,wlter;
625   PetscErrorCode ierr;
626 
627   PetscFunctionBegin;
628   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
629   /* Cannot compute LTE in first step or in restart after event */
630   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
631   /* Compute LTE using backward differences with non-constant time step */
632   {
633     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
634     PetscReal   a = 1 + h_prev/h;
635     PetscScalar scal[3]; Vec vecs[3];
636     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
637     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
638     ierr = VecCopy(X,Y);CHKERRQ(ierr);
639     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
640     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
641   }
642   if (order) *order = 2;
643   PetscFunctionReturn(0);
644 }
645 
646 static PetscErrorCode TSRollBack_Theta(TS ts)
647 {
648   TS_Theta       *th = (TS_Theta*)ts->data;
649   PetscInt       ncost;
650   PetscErrorCode ierr;
651 
652   PetscFunctionBegin;
653   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
654   if (ts->vec_costintegral && ts->costintegralfwd) {
655     ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
656   }
657   th->status = TS_STEP_INCOMPLETE;
658   if (ts->mat_sensip) {
659     ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
660   }
661   if (ts->vecs_integral_sensip) {
662     for (ncost=0;ncost<ts->numcost;ncost++) {
663       ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr);
664     }
665   }
666   PetscFunctionReturn(0);
667 }
668 
669 static PetscErrorCode TSForwardStep_Theta(TS ts)
670 {
671   TS_Theta       *th = (TS_Theta*)ts->data;
672   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
673   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
674   PetscInt       ncost,ntlm;
675   KSP            ksp;
676   Mat            J,Jp;
677   PetscReal      shift;
678   PetscScalar    *barr,*xarr;
679   PetscErrorCode ierr;
680 
681   PetscFunctionBegin;
682   ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
683 
684   for (ncost=0; ncost<ts->numcost; ncost++) {
685     if (ts->vecs_integral_sensip) {
686       ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr);
687     }
688   }
689 
690   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
691   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
692 
693   /* Build RHS */
694   if (th->endpoint) { /* 2-stage method*/
695     shift = 1./((th->Theta-1.)*th->time_step);
696     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
697     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
698     ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
699 
700     /* Add the f_p forcing terms */
701     if (ts->Jacp) {
702       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
703       ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
704       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
705       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
706     }
707   } else { /* 1-stage method */
708     shift = 0.0;
709     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
710     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
711     ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr);
712 
713     /* Add the f_p forcing terms */
714     if (ts->Jacp) {
715       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
716       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
717     }
718   }
719 
720   /* Build LHS */
721   shift = 1/(th->Theta*th->time_step);
722   if (th->endpoint) {
723     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
724   } else {
725     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
726   }
727   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
728 
729   /*
730     Evaluate the first stage of integral gradients with the 2-stage method:
731     drdu|t_n*S(t_n) + drdp|t_n
732     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
733   */
734   if (th->endpoint) { /* 2-stage method only */
735     if (ts->vecs_integral_sensip) {
736       ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr);
737       if (ts->vecs_drdp) {
738         ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
739       }
740       for (ncost=0; ncost<ts->numcost; ncost++) {
741         ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr);
742         if (ts->vecs_drdp) {
743           ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr);
744         }
745         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr);
746       }
747     }
748   }
749 
750   /* Solve the tangent linear equation for forward sensitivities to parameters */
751   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
752     KSPConvergedReason kspreason;
753     ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr);
754     ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr);
755     if (th->endpoint) {
756       ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr);
757       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
758       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr);
759       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
760       ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
761     } else {
762       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr);
763     }
764     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
765     if (kspreason < 0) {
766       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
767       ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr);
768     }
769     ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr);
770     ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr);
771   }
772 
773 
774   /*
775     Evaluate the second stage of integral gradients with the 2-stage method:
776     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
777   */
778   if (ts->vecs_integral_sensip) {
779     if (!th->endpoint) {
780       ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
781       ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr);
782       if (ts->vecs_drdp) {
783         ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
784       }
785       for (ncost=0; ncost<ts->numcost; ncost++) {
786         ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr);
787         if (ts->vecs_drdp) {
788           ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr);
789         }
790         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr);
791       }
792       ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
793     } else {
794       ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr);
795       if (ts->vecs_drdp) {
796         ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
797       }
798       for (ncost=0; ncost<ts->numcost; ncost++) {
799         ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr);
800         if (ts->vecs_drdp) {
801           ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr);
802         }
803         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr);
804       }
805     }
806   } else {
807     if (!th->endpoint) {
808       ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
809     }
810   }
811   PetscFunctionReturn(0);
812 }
813 
814 static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip)
815 {
816   TS_Theta *th = (TS_Theta*)ts->data;
817 
818   PetscFunctionBegin;
819   if (ns) *ns = 1;
820   if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip);
821   PetscFunctionReturn(0);
822 }
823 
824 /*------------------------------------------------------------*/
825 static PetscErrorCode TSReset_Theta(TS ts)
826 {
827   TS_Theta       *th = (TS_Theta*)ts->data;
828   PetscErrorCode ierr;
829 
830   PetscFunctionBegin;
831   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
832   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
833   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
834   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
835 
836   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
837   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
838 
839   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
840   PetscFunctionReturn(0);
841 }
842 
843 static PetscErrorCode TSAdjointReset_Theta(TS ts)
844 {
845   TS_Theta       *th = (TS_Theta*)ts->data;
846   PetscErrorCode ierr;
847 
848   PetscFunctionBegin;
849   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
850   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
851   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
852   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
853   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
854   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
855   PetscFunctionReturn(0);
856 }
857 
858 static PetscErrorCode TSDestroy_Theta(TS ts)
859 {
860   PetscErrorCode ierr;
861 
862   PetscFunctionBegin;
863   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
864   if (ts->dm) {
865     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
866     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
867   }
868   ierr = PetscFree(ts->data);CHKERRQ(ierr);
869   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
870   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
871   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
872   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
873   PetscFunctionReturn(0);
874 }
875 
876 /*
877   This defines the nonlinear equation that is to be solved with SNES
878   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
879 */
880 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
881 {
882   TS_Theta       *th = (TS_Theta*)ts->data;
883   PetscErrorCode ierr;
884   Vec            X0,Xdot;
885   DM             dm,dmsave;
886   PetscReal      shift = 1/(th->Theta*ts->time_step);
887 
888   PetscFunctionBegin;
889   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
890   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
891   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
892   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
893 
894   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
895   dmsave = ts->dm;
896   ts->dm = dm;
897   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
898   ts->dm = dmsave;
899   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
900   PetscFunctionReturn(0);
901 }
902 
903 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
904 {
905   TS_Theta       *th = (TS_Theta*)ts->data;
906   PetscErrorCode ierr;
907   Vec            Xdot;
908   DM             dm,dmsave;
909   PetscReal      shift = 1/(th->Theta*ts->time_step);
910 
911   PetscFunctionBegin;
912   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
913   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
914   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
915 
916   dmsave = ts->dm;
917   ts->dm = dm;
918   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
919   ts->dm = dmsave;
920   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
921   PetscFunctionReturn(0);
922 }
923 
924 static PetscErrorCode TSForwardSetUp_Theta(TS ts)
925 {
926   TS_Theta       *th = (TS_Theta*)ts->data;
927   PetscErrorCode ierr;
928 
929   PetscFunctionBegin;
930   /* combine sensitivities to parameters and sensitivities to initial values into one array */
931   th->num_tlm = ts->num_parameters;
932   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr);
933   if (ts->vecs_integral_sensip) {
934     ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr);
935   }
936   /* backup sensitivity results for roll-backs */
937   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr);
938 
939   if (ts->vecs_integral_sensip) {
940     ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
941   }
942   ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
943   PetscFunctionReturn(0);
944 }
945 
946 static PetscErrorCode TSForwardReset_Theta(TS ts)
947 {
948   TS_Theta       *th = (TS_Theta*)ts->data;
949   PetscErrorCode ierr;
950 
951   PetscFunctionBegin;
952   if (ts->vecs_integral_sensip) {
953     ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr);
954     ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
955   }
956   ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
957   ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr);
958   ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr);
959   PetscFunctionReturn(0);
960 }
961 
962 static PetscErrorCode TSSetUp_Theta(TS ts)
963 {
964   TS_Theta       *th = (TS_Theta*)ts->data;
965   PetscBool      match;
966   PetscErrorCode ierr;
967 
968   PetscFunctionBegin;
969   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
970     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
971   }
972   if (!th->X) {
973     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
974   }
975   if (!th->Xdot) {
976     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
977   }
978   if (!th->X0) {
979     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
980   }
981   if (th->endpoint) {
982     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
983   }
984 
985   th->order = (th->Theta == 0.5) ? 2 : 1;
986 
987   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
988   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
989   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
990 
991   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
992   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
993   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
994   if (!match) {
995     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
996     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
997   }
998   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
999   PetscFunctionReturn(0);
1000 }
1001 
1002 /*------------------------------------------------------------*/
1003 
1004 static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1005 {
1006   TS_Theta       *th = (TS_Theta*)ts->data;
1007   PetscErrorCode ierr;
1008 
1009   PetscFunctionBegin;
1010   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
1011   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
1012   if (ts->vecs_sensip) {
1013     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
1014   }
1015   if (ts->vecs_sensi2) {
1016     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
1017     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
1018     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
1019     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
1020     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1021   }
1022   if (ts->vecs_sensi2p) {
1023     ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
1024     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
1025     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
1026     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1027   }
1028   PetscFunctionReturn(0);
1029 }
1030 
1031 static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1032 {
1033   TS_Theta       *th = (TS_Theta*)ts->data;
1034   PetscErrorCode ierr;
1035 
1036   PetscFunctionBegin;
1037   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
1038   {
1039     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
1040     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
1041     ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
1042   }
1043   ierr = PetscOptionsTail();CHKERRQ(ierr);
1044   PetscFunctionReturn(0);
1045 }
1046 
1047 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1048 {
1049   TS_Theta       *th = (TS_Theta*)ts->data;
1050   PetscBool      iascii;
1051   PetscErrorCode ierr;
1052 
1053   PetscFunctionBegin;
1054   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1055   if (iascii) {
1056     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
1057     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
1058   }
1059   PetscFunctionReturn(0);
1060 }
1061 
1062 static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
1063 {
1064   TS_Theta *th = (TS_Theta*)ts->data;
1065 
1066   PetscFunctionBegin;
1067   *theta = th->Theta;
1068   PetscFunctionReturn(0);
1069 }
1070 
1071 static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
1072 {
1073   TS_Theta *th = (TS_Theta*)ts->data;
1074 
1075   PetscFunctionBegin;
1076   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
1077   th->Theta = theta;
1078   th->order = (th->Theta == 0.5) ? 2 : 1;
1079   PetscFunctionReturn(0);
1080 }
1081 
1082 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
1083 {
1084   TS_Theta *th = (TS_Theta*)ts->data;
1085 
1086   PetscFunctionBegin;
1087   *endpoint = th->endpoint;
1088   PetscFunctionReturn(0);
1089 }
1090 
1091 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1092 {
1093   TS_Theta *th = (TS_Theta*)ts->data;
1094 
1095   PetscFunctionBegin;
1096   th->endpoint = flg;
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 #if defined(PETSC_HAVE_COMPLEX)
1101 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1102 {
1103   PetscComplex   z   = xr + xi*PETSC_i,f;
1104   TS_Theta       *th = (TS_Theta*)ts->data;
1105   const PetscReal one = 1.0;
1106 
1107   PetscFunctionBegin;
1108   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1109   *yr = PetscRealPartComplex(f);
1110   *yi = PetscImaginaryPartComplex(f);
1111   PetscFunctionReturn(0);
1112 }
1113 #endif
1114 
1115 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
1116 {
1117   TS_Theta     *th = (TS_Theta*)ts->data;
1118 
1119   PetscFunctionBegin;
1120   if (ns) *ns = 1;
1121   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 /* ------------------------------------------------------------ */
1126 /*MC
1127       TSTHETA - DAE solver using the implicit Theta method
1128 
1129    Level: beginner
1130 
1131    Options Database:
1132 +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1133 .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
1134 -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
1135 
1136    Notes:
1137 $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
1138 $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
1139 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
1140 
1141    This method can be applied to DAE.
1142 
1143    This method is cast as a 1-stage implicit Runge-Kutta method.
1144 
1145 .vb
1146   Theta | Theta
1147   -------------
1148         |  1
1149 .ve
1150 
1151    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1152 
1153    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1154 
1155 .vb
1156   0 | 0         0
1157   1 | 1-Theta   Theta
1158   -------------------
1159     | 1-Theta   Theta
1160 .ve
1161 
1162    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1163 
1164    To apply a diagonally implicit RK method to DAE, the stage formula
1165 
1166 $  Y_i = X + h sum_j a_ij Y'_j
1167 
1168    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1169 
1170 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1171 
1172 M*/
1173 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1174 {
1175   TS_Theta       *th;
1176   PetscErrorCode ierr;
1177 
1178   PetscFunctionBegin;
1179   ts->ops->reset           = TSReset_Theta;
1180   ts->ops->adjointreset    = TSAdjointReset_Theta;
1181   ts->ops->destroy         = TSDestroy_Theta;
1182   ts->ops->view            = TSView_Theta;
1183   ts->ops->setup           = TSSetUp_Theta;
1184   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1185   ts->ops->adjointreset    = TSAdjointReset_Theta;
1186   ts->ops->step            = TSStep_Theta;
1187   ts->ops->interpolate     = TSInterpolate_Theta;
1188   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
1189   ts->ops->rollback        = TSRollBack_Theta;
1190   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
1191   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
1192   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1193 #if defined(PETSC_HAVE_COMPLEX)
1194   ts->ops->linearstability = TSComputeLinearStability_Theta;
1195 #endif
1196   ts->ops->getstages       = TSGetStages_Theta;
1197   ts->ops->adjointstep     = TSAdjointStep_Theta;
1198   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1199   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
1200   ts->default_adapt_type   = TSADAPTNONE;
1201 
1202   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
1203   ts->ops->forwardreset     = TSForwardReset_Theta;
1204   ts->ops->forwardstep      = TSForwardStep_Theta;
1205   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1206 
1207   ts->usessnes = PETSC_TRUE;
1208 
1209   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1210   ts->data = (void*)th;
1211 
1212   th->VecsDeltaLam    = NULL;
1213   th->VecsDeltaMu     = NULL;
1214   th->VecsSensiTemp   = NULL;
1215   th->VecsSensi2Temp  = NULL;
1216 
1217   th->extrapolate = PETSC_FALSE;
1218   th->Theta       = 0.5;
1219   th->order       = 2;
1220   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1221   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1222   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1223   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 /*@
1228   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
1229 
1230   Not Collective
1231 
1232   Input Parameter:
1233 .  ts - timestepping context
1234 
1235   Output Parameter:
1236 .  theta - stage abscissa
1237 
1238   Note:
1239   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
1240 
1241   Level: Advanced
1242 
1243 .seealso: TSThetaSetTheta()
1244 @*/
1245 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
1246 {
1247   PetscErrorCode ierr;
1248 
1249   PetscFunctionBegin;
1250   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1251   PetscValidPointer(theta,2);
1252   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
1253   PetscFunctionReturn(0);
1254 }
1255 
1256 /*@
1257   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
1258 
1259   Not Collective
1260 
1261   Input Parameter:
1262 +  ts - timestepping context
1263 -  theta - stage abscissa
1264 
1265   Options Database:
1266 .  -ts_theta_theta <theta>
1267 
1268   Level: Intermediate
1269 
1270 .seealso: TSThetaGetTheta()
1271 @*/
1272 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
1273 {
1274   PetscErrorCode ierr;
1275 
1276   PetscFunctionBegin;
1277   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1278   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
1279   PetscFunctionReturn(0);
1280 }
1281 
1282 /*@
1283   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1284 
1285   Not Collective
1286 
1287   Input Parameter:
1288 .  ts - timestepping context
1289 
1290   Output Parameter:
1291 .  endpoint - PETSC_TRUE when using the endpoint variant
1292 
1293   Level: Advanced
1294 
1295 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
1296 @*/
1297 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
1298 {
1299   PetscErrorCode ierr;
1300 
1301   PetscFunctionBegin;
1302   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1303   PetscValidPointer(endpoint,2);
1304   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
1305   PetscFunctionReturn(0);
1306 }
1307 
1308 /*@
1309   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1310 
1311   Not Collective
1312 
1313   Input Parameter:
1314 +  ts - timestepping context
1315 -  flg - PETSC_TRUE to use the endpoint variant
1316 
1317   Options Database:
1318 .  -ts_theta_endpoint <flg>
1319 
1320   Level: Intermediate
1321 
1322 .seealso: TSTHETA, TSCN
1323 @*/
1324 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1325 {
1326   PetscErrorCode ierr;
1327 
1328   PetscFunctionBegin;
1329   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1330   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1331   PetscFunctionReturn(0);
1332 }
1333 
1334 /*
1335  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1336  * The creation functions for these specializations are below.
1337  */
1338 
1339 static PetscErrorCode TSSetUp_BEuler(TS ts)
1340 {
1341   TS_Theta       *th = (TS_Theta*)ts->data;
1342   PetscErrorCode ierr;
1343 
1344   PetscFunctionBegin;
1345   if (th->Theta != 1.0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (1) of theta when using backward Euler\n");
1346   if (th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the endpoint form of the Theta methods when using backward Euler\n");
1347   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1348   PetscFunctionReturn(0);
1349 }
1350 
1351 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1352 {
1353   PetscFunctionBegin;
1354   PetscFunctionReturn(0);
1355 }
1356 
1357 /*MC
1358       TSBEULER - ODE solver using the implicit backward Euler method
1359 
1360   Level: beginner
1361 
1362   Notes:
1363   TSBEULER is equivalent to TSTHETA with Theta=1.0
1364 
1365 $  -ts_type theta -ts_theta_theta 1.0
1366 
1367 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1368 
1369 M*/
1370 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1371 {
1372   PetscErrorCode ierr;
1373 
1374   PetscFunctionBegin;
1375   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1376   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
1377   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
1378   ts->ops->setup = TSSetUp_BEuler;
1379   ts->ops->view  = TSView_BEuler;
1380   PetscFunctionReturn(0);
1381 }
1382 
1383 static PetscErrorCode TSSetUp_CN(TS ts)
1384 {
1385   TS_Theta       *th = (TS_Theta*)ts->data;
1386   PetscErrorCode ierr;
1387 
1388   PetscFunctionBegin;
1389   if (th->Theta != 0.5) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (0.5) of theta when using Crank-Nicolson\n");
1390   if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the midpoint form of the Theta methods when using Crank-Nicolson\n");
1391   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1392   PetscFunctionReturn(0);
1393 }
1394 
1395 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1396 {
1397   PetscFunctionBegin;
1398   PetscFunctionReturn(0);
1399 }
1400 
1401 /*MC
1402       TSCN - ODE solver using the implicit Crank-Nicolson method.
1403 
1404   Level: beginner
1405 
1406   Notes:
1407   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
1408 
1409 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1410 
1411 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1412 
1413 M*/
1414 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1415 {
1416   PetscErrorCode ierr;
1417 
1418   PetscFunctionBegin;
1419   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1420   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1421   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
1422   ts->ops->setup = TSSetUp_CN;
1423   ts->ops->view  = TSView_CN;
1424   PetscFunctionReturn(0);
1425 }
1426