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