xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision b2d3fc2fff6b0b3e3ee1d03a12efee565de8bc1c)
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    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 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
945 {
946   TS_Theta       *th = (TS_Theta*)ts->data;
947   PetscErrorCode ierr;
948   Vec            X0,Xdot;
949   DM             dm,dmsave;
950   PetscReal      shift = th->shift;
951 
952   PetscFunctionBegin;
953   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
954   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
955   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
956   if (x != X0) {
957     ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
958   } else {
959     ierr = VecZeroEntries(Xdot);CHKERRQ(ierr);
960   }
961   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
962   dmsave = ts->dm;
963   ts->dm = dm;
964   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
965   ts->dm = dmsave;
966   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
967   PetscFunctionReturn(0);
968 }
969 
970 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
971 {
972   TS_Theta       *th = (TS_Theta*)ts->data;
973   PetscErrorCode ierr;
974   Vec            Xdot;
975   DM             dm,dmsave;
976   PetscReal      shift = th->shift;
977 
978   PetscFunctionBegin;
979   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
980   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
981   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
982 
983   dmsave = ts->dm;
984   ts->dm = dm;
985   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
986   ts->dm = dmsave;
987   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
988   PetscFunctionReturn(0);
989 }
990 
991 static PetscErrorCode TSForwardSetUp_Theta(TS ts)
992 {
993   TS_Theta       *th = (TS_Theta*)ts->data;
994   TS             quadts = ts->quadraturets;
995   PetscErrorCode ierr;
996 
997   PetscFunctionBegin;
998   /* combine sensitivities to parameters and sensitivities to initial values into one array */
999   th->num_tlm = ts->num_parameters;
1000   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr);
1001   if (quadts && quadts->mat_sensip) {
1002     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
1003     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr);
1004   }
1005   /* backup sensitivity results for roll-backs */
1006   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr);
1007 
1008   ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 static PetscErrorCode TSForwardReset_Theta(TS ts)
1013 {
1014   TS_Theta       *th = (TS_Theta*)ts->data;
1015   TS             quadts = ts->quadraturets;
1016   PetscErrorCode ierr;
1017 
1018   PetscFunctionBegin;
1019   if (quadts && quadts->mat_sensip) {
1020     ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr);
1021     ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr);
1022   }
1023   ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
1024   ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr);
1025   ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr);
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 static PetscErrorCode TSSetUp_Theta(TS ts)
1030 {
1031   TS_Theta       *th = (TS_Theta*)ts->data;
1032   TS             quadts = ts->quadraturets;
1033   PetscBool      match;
1034   PetscErrorCode ierr;
1035 
1036   PetscFunctionBegin;
1037   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
1038     ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr);
1039   }
1040   if (!th->X) {
1041     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
1042   }
1043   if (!th->Xdot) {
1044     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
1045   }
1046   if (!th->X0) {
1047     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
1048   }
1049   if (th->endpoint) {
1050     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
1051   }
1052 
1053   th->order = (th->Theta == 0.5) ? 2 : 1;
1054   th->shift = 1/(th->Theta*ts->time_step);
1055 
1056   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
1057   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
1058   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
1059 
1060   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
1061   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
1062   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
1063   if (!match) {
1064     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
1065     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
1066   }
1067   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
1068   PetscFunctionReturn(0);
1069 }
1070 
1071 /*------------------------------------------------------------*/
1072 
1073 static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1074 {
1075   TS_Theta       *th = (TS_Theta*)ts->data;
1076   PetscErrorCode ierr;
1077 
1078   PetscFunctionBegin;
1079   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
1080   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
1081   if (ts->vecs_sensip) {
1082     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
1083   }
1084   if (ts->vecs_sensi2) {
1085     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
1086     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
1087     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
1088     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
1089     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1090   }
1091   if (ts->vecs_sensi2p) {
1092     ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
1093     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
1094     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
1095     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1096   }
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1101 {
1102   TS_Theta       *th = (TS_Theta*)ts->data;
1103   PetscErrorCode ierr;
1104 
1105   PetscFunctionBegin;
1106   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
1107   {
1108     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
1109     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
1110     ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
1111   }
1112   ierr = PetscOptionsTail();CHKERRQ(ierr);
1113   PetscFunctionReturn(0);
1114 }
1115 
1116 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1117 {
1118   TS_Theta       *th = (TS_Theta*)ts->data;
1119   PetscBool      iascii;
1120   PetscErrorCode ierr;
1121 
1122   PetscFunctionBegin;
1123   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1124   if (iascii) {
1125     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
1126     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
1127   }
1128   PetscFunctionReturn(0);
1129 }
1130 
1131 static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
1132 {
1133   TS_Theta *th = (TS_Theta*)ts->data;
1134 
1135   PetscFunctionBegin;
1136   *theta = th->Theta;
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
1141 {
1142   TS_Theta *th = (TS_Theta*)ts->data;
1143 
1144   PetscFunctionBegin;
1145   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
1146   th->Theta = theta;
1147   th->order = (th->Theta == 0.5) ? 2 : 1;
1148   PetscFunctionReturn(0);
1149 }
1150 
1151 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
1152 {
1153   TS_Theta *th = (TS_Theta*)ts->data;
1154 
1155   PetscFunctionBegin;
1156   *endpoint = th->endpoint;
1157   PetscFunctionReturn(0);
1158 }
1159 
1160 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1161 {
1162   TS_Theta *th = (TS_Theta*)ts->data;
1163 
1164   PetscFunctionBegin;
1165   th->endpoint = flg;
1166   PetscFunctionReturn(0);
1167 }
1168 
1169 #if defined(PETSC_HAVE_COMPLEX)
1170 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1171 {
1172   PetscComplex   z   = xr + xi*PETSC_i,f;
1173   TS_Theta       *th = (TS_Theta*)ts->data;
1174   const PetscReal one = 1.0;
1175 
1176   PetscFunctionBegin;
1177   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1178   *yr = PetscRealPartComplex(f);
1179   *yi = PetscImaginaryPartComplex(f);
1180   PetscFunctionReturn(0);
1181 }
1182 #endif
1183 
1184 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
1185 {
1186   TS_Theta     *th = (TS_Theta*)ts->data;
1187 
1188   PetscFunctionBegin;
1189   if (ns) *ns = 1;
1190   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
1191   PetscFunctionReturn(0);
1192 }
1193 
1194 /* ------------------------------------------------------------ */
1195 /*MC
1196       TSTHETA - DAE solver using the implicit Theta method
1197 
1198    Level: beginner
1199 
1200    Options Database:
1201 +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1202 .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
1203 -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
1204 
1205    Notes:
1206 $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
1207 $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
1208 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
1209 
1210    This method can be applied to DAE.
1211 
1212    This method is cast as a 1-stage implicit Runge-Kutta method.
1213 
1214 .vb
1215   Theta | Theta
1216   -------------
1217         |  1
1218 .ve
1219 
1220    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1221 
1222    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1223 
1224 .vb
1225   0 | 0         0
1226   1 | 1-Theta   Theta
1227   -------------------
1228     | 1-Theta   Theta
1229 .ve
1230 
1231    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1232 
1233    To apply a diagonally implicit RK method to DAE, the stage formula
1234 
1235 $  Y_i = X + h sum_j a_ij Y'_j
1236 
1237    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1238 
1239 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1240 
1241 M*/
1242 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1243 {
1244   TS_Theta       *th;
1245   PetscErrorCode ierr;
1246 
1247   PetscFunctionBegin;
1248   ts->ops->reset           = TSReset_Theta;
1249   ts->ops->adjointreset    = TSAdjointReset_Theta;
1250   ts->ops->destroy         = TSDestroy_Theta;
1251   ts->ops->view            = TSView_Theta;
1252   ts->ops->setup           = TSSetUp_Theta;
1253   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1254   ts->ops->adjointreset    = TSAdjointReset_Theta;
1255   ts->ops->step            = TSStep_Theta;
1256   ts->ops->interpolate     = TSInterpolate_Theta;
1257   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
1258   ts->ops->rollback        = TSRollBack_Theta;
1259   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
1260   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
1261   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1262 #if defined(PETSC_HAVE_COMPLEX)
1263   ts->ops->linearstability = TSComputeLinearStability_Theta;
1264 #endif
1265   ts->ops->getstages       = TSGetStages_Theta;
1266   ts->ops->adjointstep     = TSAdjointStep_Theta;
1267   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1268   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
1269   ts->default_adapt_type   = TSADAPTNONE;
1270 
1271   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
1272   ts->ops->forwardreset     = TSForwardReset_Theta;
1273   ts->ops->forwardstep      = TSForwardStep_Theta;
1274   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1275 
1276   ts->usessnes = PETSC_TRUE;
1277 
1278   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1279   ts->data = (void*)th;
1280 
1281   th->VecsDeltaLam    = NULL;
1282   th->VecsDeltaMu     = NULL;
1283   th->VecsSensiTemp   = NULL;
1284   th->VecsSensi2Temp  = NULL;
1285 
1286   th->extrapolate = PETSC_FALSE;
1287   th->Theta       = 0.5;
1288   th->order       = 2;
1289   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1290   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1291   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1292   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1293   PetscFunctionReturn(0);
1294 }
1295 
1296 /*@
1297   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
1298 
1299   Not Collective
1300 
1301   Input Parameter:
1302 .  ts - timestepping context
1303 
1304   Output Parameter:
1305 .  theta - stage abscissa
1306 
1307   Note:
1308   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
1309 
1310   Level: Advanced
1311 
1312 .seealso: TSThetaSetTheta()
1313 @*/
1314 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
1315 {
1316   PetscErrorCode ierr;
1317 
1318   PetscFunctionBegin;
1319   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1320   PetscValidPointer(theta,2);
1321   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
1322   PetscFunctionReturn(0);
1323 }
1324 
1325 /*@
1326   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
1327 
1328   Not Collective
1329 
1330   Input Parameter:
1331 +  ts - timestepping context
1332 -  theta - stage abscissa
1333 
1334   Options Database:
1335 .  -ts_theta_theta <theta>
1336 
1337   Level: Intermediate
1338 
1339 .seealso: TSThetaGetTheta()
1340 @*/
1341 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
1342 {
1343   PetscErrorCode ierr;
1344 
1345   PetscFunctionBegin;
1346   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1347   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
1348   PetscFunctionReturn(0);
1349 }
1350 
1351 /*@
1352   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1353 
1354   Not Collective
1355 
1356   Input Parameter:
1357 .  ts - timestepping context
1358 
1359   Output Parameter:
1360 .  endpoint - PETSC_TRUE when using the endpoint variant
1361 
1362   Level: Advanced
1363 
1364 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
1365 @*/
1366 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
1367 {
1368   PetscErrorCode ierr;
1369 
1370   PetscFunctionBegin;
1371   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1372   PetscValidPointer(endpoint,2);
1373   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
1374   PetscFunctionReturn(0);
1375 }
1376 
1377 /*@
1378   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1379 
1380   Not Collective
1381 
1382   Input Parameter:
1383 +  ts - timestepping context
1384 -  flg - PETSC_TRUE to use the endpoint variant
1385 
1386   Options Database:
1387 .  -ts_theta_endpoint <flg>
1388 
1389   Level: Intermediate
1390 
1391 .seealso: TSTHETA, TSCN
1392 @*/
1393 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1394 {
1395   PetscErrorCode ierr;
1396 
1397   PetscFunctionBegin;
1398   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1399   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1400   PetscFunctionReturn(0);
1401 }
1402 
1403 /*
1404  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1405  * The creation functions for these specializations are below.
1406  */
1407 
1408 static PetscErrorCode TSSetUp_BEuler(TS ts)
1409 {
1410   TS_Theta       *th = (TS_Theta*)ts->data;
1411   PetscErrorCode ierr;
1412 
1413   PetscFunctionBegin;
1414   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");
1415   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");
1416   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1417   PetscFunctionReturn(0);
1418 }
1419 
1420 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1421 {
1422   PetscFunctionBegin;
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 /*MC
1427       TSBEULER - ODE solver using the implicit backward Euler method
1428 
1429   Level: beginner
1430 
1431   Notes:
1432   TSBEULER is equivalent to TSTHETA with Theta=1.0
1433 
1434 $  -ts_type theta -ts_theta_theta 1.0
1435 
1436 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1437 
1438 M*/
1439 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1440 {
1441   PetscErrorCode ierr;
1442 
1443   PetscFunctionBegin;
1444   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1445   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
1446   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
1447   ts->ops->setup = TSSetUp_BEuler;
1448   ts->ops->view  = TSView_BEuler;
1449   PetscFunctionReturn(0);
1450 }
1451 
1452 static PetscErrorCode TSSetUp_CN(TS ts)
1453 {
1454   TS_Theta       *th = (TS_Theta*)ts->data;
1455   PetscErrorCode ierr;
1456 
1457   PetscFunctionBegin;
1458   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");
1459   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");
1460   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1461   PetscFunctionReturn(0);
1462 }
1463 
1464 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1465 {
1466   PetscFunctionBegin;
1467   PetscFunctionReturn(0);
1468 }
1469 
1470 /*MC
1471       TSCN - ODE solver using the implicit Crank-Nicolson method.
1472 
1473   Level: beginner
1474 
1475   Notes:
1476   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
1477 
1478 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1479 
1480 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1481 
1482 M*/
1483 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1484 {
1485   PetscErrorCode ierr;
1486 
1487   PetscFunctionBegin;
1488   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1489   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1490   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
1491   ts->ops->setup = TSSetUp_CN;
1492   ts->ops->view  = TSView_CN;
1493   PetscFunctionReturn(0);
1494 }
1495