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