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