xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision fecfb7146cb415f81e8cc7ef8d4ae968f3ac4b11)
1316643e7SJed Brown  /*
2316643e7SJed Brown    Code for timestepping with implicit Theta method
3316643e7SJed Brown  */
4af0996ceSBarry Smith  #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
5a055c8caSBarry Smith  #include <petscsnes.h>
61e25c274SJed Brown  #include <petscdm.h>
73c54f92cSHong Zhang  #include <petscmat.h>
8316643e7SJed Brown 
9316643e7SJed Brown  typedef struct {
101566a47fSLisandro Dalcin    PetscReal    stage_time;
111566a47fSLisandro Dalcin    Vec          X0,X,Xdot;                /* Storage for stages and time derivative */
121566a47fSLisandro Dalcin    Vec          affine;                   /* Affine vector needed for residual at beginning of step in endpoint formulation */
131566a47fSLisandro Dalcin 
141566a47fSLisandro Dalcin    PetscReal    Theta;
151566a47fSLisandro Dalcin    PetscInt     order;
161566a47fSLisandro Dalcin    PetscBool    endpoint;
171566a47fSLisandro Dalcin    PetscBool    extrapolate;
181566a47fSLisandro Dalcin 
19b5e0532dSHong Zhang    Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
20b5e0532dSHong Zhang    Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
21b5e0532dSHong Zhang    Vec          *VecsSensiTemp;           /* Vector to be timed with Jacobian transpose */
221566a47fSLisandro Dalcin    Vec          VecCostIntegral0;         /* Backup for roll-backs due to events */
23b5e0532dSHong Zhang    PetscReal    ptime;
24b1cb36f3SHong Zhang    PetscReal    time_step;
251566a47fSLisandro Dalcin 
261566a47fSLisandro Dalcin    PetscBool    adapt;                    /* Use time-step adaptivity ? */
271566a47fSLisandro Dalcin    Vec          vec_sol_prev;
281566a47fSLisandro Dalcin    Vec          vec_lte_work;
291566a47fSLisandro Dalcin 
301566a47fSLisandro Dalcin    TSStepStatus status;
31316643e7SJed Brown  } TS_Theta;
32316643e7SJed Brown 
33316643e7SJed Brown  #undef __FUNCT__
347445fe48SJed Brown  #define __FUNCT__ "TSThetaGetX0AndXdot"
357445fe48SJed Brown  static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
367445fe48SJed Brown  {
377445fe48SJed Brown    TS_Theta       *th = (TS_Theta*)ts->data;
387445fe48SJed Brown    PetscErrorCode ierr;
397445fe48SJed Brown 
407445fe48SJed Brown    PetscFunctionBegin;
417445fe48SJed Brown    if (X0) {
427445fe48SJed Brown      if (dm && dm != ts->dm) {
430d0b770aSPeter Brune        ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
447445fe48SJed Brown      } else *X0 = ts->vec_sol;
457445fe48SJed Brown    }
467445fe48SJed Brown    if (Xdot) {
477445fe48SJed Brown      if (dm && dm != ts->dm) {
480d0b770aSPeter Brune        ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
497445fe48SJed Brown      } else *Xdot = th->Xdot;
507445fe48SJed Brown    }
517445fe48SJed Brown    PetscFunctionReturn(0);
527445fe48SJed Brown  }
537445fe48SJed Brown 
540d0b770aSPeter Brune  #undef __FUNCT__
550d0b770aSPeter Brune  #define __FUNCT__ "TSThetaRestoreX0AndXdot"
560d0b770aSPeter Brune  static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
570d0b770aSPeter Brune  {
580d0b770aSPeter Brune    PetscErrorCode ierr;
590d0b770aSPeter Brune 
600d0b770aSPeter Brune    PetscFunctionBegin;
610d0b770aSPeter Brune    if (X0) {
620d0b770aSPeter Brune      if (dm && dm != ts->dm) {
630d0b770aSPeter Brune        ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
640d0b770aSPeter Brune      }
650d0b770aSPeter Brune    }
660d0b770aSPeter Brune    if (Xdot) {
670d0b770aSPeter Brune      if (dm && dm != ts->dm) {
680d0b770aSPeter Brune        ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
690d0b770aSPeter Brune      }
700d0b770aSPeter Brune    }
710d0b770aSPeter Brune    PetscFunctionReturn(0);
720d0b770aSPeter Brune  }
730d0b770aSPeter Brune 
747445fe48SJed Brown  #undef __FUNCT__
757445fe48SJed Brown  #define __FUNCT__ "DMCoarsenHook_TSTheta"
767445fe48SJed Brown  static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
777445fe48SJed Brown  {
787445fe48SJed Brown 
797445fe48SJed Brown    PetscFunctionBegin;
807445fe48SJed Brown    PetscFunctionReturn(0);
817445fe48SJed Brown  }
827445fe48SJed Brown 
837445fe48SJed Brown  #undef __FUNCT__
847445fe48SJed Brown  #define __FUNCT__ "DMRestrictHook_TSTheta"
857445fe48SJed Brown  static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
867445fe48SJed Brown  {
877445fe48SJed Brown    TS             ts = (TS)ctx;
887445fe48SJed Brown    PetscErrorCode ierr;
897445fe48SJed Brown    Vec            X0,Xdot,X0_c,Xdot_c;
907445fe48SJed Brown 
917445fe48SJed Brown    PetscFunctionBegin;
927445fe48SJed Brown    ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
937445fe48SJed Brown    ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
947445fe48SJed Brown    ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
957445fe48SJed Brown    ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
967445fe48SJed Brown    ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
977445fe48SJed Brown    ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
980d0b770aSPeter Brune    ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
990d0b770aSPeter Brune    ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
1007445fe48SJed Brown    PetscFunctionReturn(0);
1017445fe48SJed Brown  }
1027445fe48SJed Brown 
1037445fe48SJed Brown  #undef __FUNCT__
104258e1594SPeter Brune  #define __FUNCT__ "DMSubDomainHook_TSTheta"
105258e1594SPeter Brune  static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
106258e1594SPeter Brune  {
107258e1594SPeter Brune 
108258e1594SPeter Brune    PetscFunctionBegin;
109258e1594SPeter Brune    PetscFunctionReturn(0);
110258e1594SPeter Brune  }
111258e1594SPeter Brune 
112258e1594SPeter Brune  #undef __FUNCT__
113258e1594SPeter Brune  #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta"
114258e1594SPeter Brune  static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
115258e1594SPeter Brune  {
116258e1594SPeter Brune    TS             ts = (TS)ctx;
117258e1594SPeter Brune    PetscErrorCode ierr;
118258e1594SPeter Brune    Vec            X0,Xdot,X0_sub,Xdot_sub;
119258e1594SPeter Brune 
120258e1594SPeter Brune    PetscFunctionBegin;
121258e1594SPeter Brune    ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
122258e1594SPeter Brune    ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
123258e1594SPeter Brune 
124258e1594SPeter Brune    ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
125258e1594SPeter Brune    ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126258e1594SPeter Brune 
127258e1594SPeter Brune    ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
128258e1594SPeter Brune    ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
129258e1594SPeter Brune 
130258e1594SPeter Brune    ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
131258e1594SPeter Brune    ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
132258e1594SPeter Brune    PetscFunctionReturn(0);
133258e1594SPeter Brune  }
134258e1594SPeter Brune 
1353b1890cdSShri Abhyankar  #undef __FUNCT__
136b1cb36f3SHong Zhang  #define __FUNCT__ "TSForwardCostIntegral_Theta"
137b1cb36f3SHong Zhang  static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
138b1cb36f3SHong Zhang  {
139b1cb36f3SHong Zhang    TS_Theta       *th = (TS_Theta*)ts->data;
140b1cb36f3SHong Zhang    PetscErrorCode ierr;
141b1cb36f3SHong Zhang 
142b1cb36f3SHong Zhang    PetscFunctionBegin;
143b1cb36f3SHong Zhang    /* backup cost integral */
144b1cb36f3SHong Zhang    ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr);
145b1cb36f3SHong Zhang    if (th->endpoint) {
146b1cb36f3SHong Zhang      ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
1471566a47fSLisandro Dalcin      ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
148b1cb36f3SHong Zhang    }
149b1cb36f3SHong Zhang    ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
150b1cb36f3SHong Zhang    if (th->endpoint) {
151b1cb36f3SHong Zhang      ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
152b1cb36f3SHong Zhang    } else {
153b1cb36f3SHong Zhang      ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
154b1cb36f3SHong Zhang    }
155b1cb36f3SHong Zhang    PetscFunctionReturn(0);
156b1cb36f3SHong Zhang  }
157b1cb36f3SHong Zhang 
158b1cb36f3SHong Zhang  #undef __FUNCT__
159b1cb36f3SHong Zhang  #define __FUNCT__ "TSAdjointCostIntegral_Theta"
160b1cb36f3SHong Zhang  static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
161b1cb36f3SHong Zhang  {
162b1cb36f3SHong Zhang    TS_Theta       *th = (TS_Theta*)ts->data;
163b1cb36f3SHong Zhang    PetscErrorCode ierr;
164b1cb36f3SHong Zhang 
165b1cb36f3SHong Zhang    PetscFunctionBegin;
166b1cb36f3SHong Zhang    if (th->endpoint) {
167b1cb36f3SHong Zhang      /* Evolve ts->vec_costintegral to compute integrals */
168b1cb36f3SHong Zhang      ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
169b1cb36f3SHong Zhang      ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
170b1cb36f3SHong Zhang      if (th->Theta!=1) {
171b1cb36f3SHong Zhang        ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
1721566a47fSLisandro Dalcin        ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1),ts->vec_costintegrand);CHKERRQ(ierr);
173b1cb36f3SHong Zhang      }
174b1cb36f3SHong Zhang    }else {
175b1cb36f3SHong Zhang      ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
176b1cb36f3SHong Zhang      ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
177b1cb36f3SHong Zhang    }
17824655328SShri    PetscFunctionReturn(0);
17924655328SShri  }
18024655328SShri 
18124655328SShri  #undef __FUNCT__
1821566a47fSLisandro Dalcin  #define __FUNCT__ "TS_SNESSolve"
1831566a47fSLisandro Dalcin  static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
1841566a47fSLisandro Dalcin  {
1851566a47fSLisandro Dalcin    PetscInt       nits,lits;
1861566a47fSLisandro Dalcin    PetscErrorCode ierr;
1871566a47fSLisandro Dalcin 
1881566a47fSLisandro Dalcin    PetscFunctionBegin;
1891566a47fSLisandro Dalcin    ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
1901566a47fSLisandro Dalcin    ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
1911566a47fSLisandro Dalcin    ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1921566a47fSLisandro Dalcin    ts->snes_its += nits; ts->ksp_its += lits;
1931566a47fSLisandro Dalcin    PetscFunctionReturn(0);
1941566a47fSLisandro Dalcin  }
1951566a47fSLisandro Dalcin 
1961566a47fSLisandro Dalcin  #undef __FUNCT__
197316643e7SJed Brown  #define __FUNCT__ "TSStep_Theta"
198193ac0bcSJed Brown  static PetscErrorCode TSStep_Theta(TS ts)
199316643e7SJed Brown  {
200316643e7SJed Brown    TS_Theta       *th = (TS_Theta*)ts->data;
2011566a47fSLisandro Dalcin    PetscInt       rejections = 0;
2024957b756SLisandro Dalcin    PetscBool      stageok,accept = PETSC_TRUE;
2031566a47fSLisandro Dalcin    PetscReal      next_time_step = ts->time_step;
204051f2191SLisandro Dalcin    PetscErrorCode ierr;
205316643e7SJed Brown 
206316643e7SJed Brown    PetscFunctionBegin;
2071566a47fSLisandro Dalcin    if (!ts->steprollback) {
2081566a47fSLisandro Dalcin      if (th->adapt) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
2093b1890cdSShri Abhyankar      ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
2101566a47fSLisandro Dalcin    }
2111566a47fSLisandro Dalcin 
2121566a47fSLisandro Dalcin    th->status = TS_STEP_INCOMPLETE;
2131566a47fSLisandro Dalcin    while (!ts->reason && th->status != TS_STEP_COMPLETE) {
2141566a47fSLisandro Dalcin 
2151566a47fSLisandro Dalcin      PetscReal shift = 1/(th->Theta*ts->time_step);
2161566a47fSLisandro Dalcin      th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
217316643e7SJed Brown 
218be5899b3SLisandro Dalcin      ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
219*fecfb714SLisandro Dalcin      if (th->extrapolate && !ts->steprestart) {
220be5899b3SLisandro Dalcin        ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr);
221be5899b3SLisandro Dalcin      }
222eb284becSJed Brown      if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
223eb284becSJed Brown        if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
2241566a47fSLisandro Dalcin        ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
2251566a47fSLisandro Dalcin        ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
2261566a47fSLisandro Dalcin        ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
2271566a47fSLisandro Dalcin      } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2281566a47fSLisandro Dalcin        ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
229eb284becSJed Brown      }
230be5899b3SLisandro Dalcin      ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
2311566a47fSLisandro Dalcin      ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
232be5899b3SLisandro Dalcin      ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
2331566a47fSLisandro Dalcin      ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
234*fecfb714SLisandro Dalcin      if (!stageok) goto reject_step;
235051f2191SLisandro Dalcin 
236051f2191SLisandro Dalcin      th->status = TS_STEP_PENDING;
2371566a47fSLisandro Dalcin      if (th->endpoint) {
2381566a47fSLisandro Dalcin        ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
2391566a47fSLisandro Dalcin      } else {
240cb3a5882SLisandro Dalcin        ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr);
2411566a47fSLisandro Dalcin        ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
2421566a47fSLisandro Dalcin      }
2431566a47fSLisandro Dalcin      ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2441566a47fSLisandro Dalcin      th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2451566a47fSLisandro Dalcin      if (!accept) {
2461566a47fSLisandro Dalcin        ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
247be5899b3SLisandro Dalcin        ts->time_step = next_time_step;
248be5899b3SLisandro Dalcin        goto reject_step;
249051f2191SLisandro Dalcin      }
2503b1890cdSShri Abhyankar 
251b1cb36f3SHong Zhang      if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
252b1cb36f3SHong Zhang        th->ptime     = ts->ptime;
253b1cb36f3SHong Zhang        th->time_step = ts->time_step;
25417145e75SHong Zhang      }
2551566a47fSLisandro Dalcin 
2562b5a38e1SLisandro Dalcin      ts->ptime += ts->time_step;
257cdbf8f93SLisandro Dalcin      ts->time_step = next_time_step;
258051f2191SLisandro Dalcin      break;
259051f2191SLisandro Dalcin 
260051f2191SLisandro Dalcin    reject_step:
261*fecfb714SLisandro Dalcin      ts->reject++; accept = PETSC_FALSE;
2621566a47fSLisandro Dalcin      if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
263051f2191SLisandro Dalcin        ts->reason = TS_DIVERGED_STEP_REJECTED;
2641566a47fSLisandro Dalcin        ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
2653b1890cdSShri Abhyankar      }
2663b1890cdSShri Abhyankar    }
267316643e7SJed Brown    PetscFunctionReturn(0);
268316643e7SJed Brown  }
269316643e7SJed Brown 
270cd652676SJed Brown  #undef __FUNCT__
27142f2b339SBarry Smith  #define __FUNCT__ "TSAdjointStep_Theta"
27242f2b339SBarry Smith  static PetscErrorCode TSAdjointStep_Theta(TS ts)
2732ca6e920SHong Zhang  {
2742ca6e920SHong Zhang    TS_Theta            *th = (TS_Theta*)ts->data;
275b5e0532dSHong Zhang    Vec                 *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
2762ca6e920SHong Zhang    PetscInt            nadj;
2772ca6e920SHong Zhang    PetscErrorCode      ierr;
2782ca6e920SHong Zhang    Mat                 J,Jp;
2792ca6e920SHong Zhang    KSP                 ksp;
2802ca6e920SHong Zhang    PetscReal           shift;
2812ca6e920SHong Zhang 
2822ca6e920SHong Zhang    PetscFunctionBegin;
2832ca6e920SHong Zhang 
2842ca6e920SHong Zhang    th->status = TS_STEP_INCOMPLETE;
285302440fdSBarry Smith    ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
2862ca6e920SHong Zhang    ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
287b5e0532dSHong Zhang 
288b5e0532dSHong Zhang    /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
2893fd52205SHong Zhang    th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
290b5e0532dSHong Zhang    th->ptime      = ts->ptime + ts->time_step;
2912ca6e920SHong Zhang 
292a4cab896SHong Zhang    /* Build RHS */
2932c39e106SBarry Smith    if (ts->vec_costintegral) { /* Cost function has an integral term */
29405755b9cSHong Zhang      if (th->endpoint) {
295d4aa0a58SBarry Smith        ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr);
29605755b9cSHong Zhang      }else {
297d4aa0a58SBarry Smith        ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
29805755b9cSHong Zhang      }
29905755b9cSHong Zhang    }
300abc2977eSBarry Smith    for (nadj=0; nadj<ts->numcost; nadj++) {
301b5e0532dSHong Zhang      ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
302b5e0532dSHong Zhang      ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
3032c39e106SBarry Smith      if (ts->vec_costintegral) {
304b5e0532dSHong Zhang        ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
30536eaed60SHong Zhang      }
3062ca6e920SHong Zhang    }
3073c54f92cSHong Zhang 
3082ca6e920SHong Zhang    /* Build LHS */
3092ca6e920SHong Zhang    shift = -1./(th->Theta*ts->time_step);
3103c54f92cSHong Zhang    if (th->endpoint) {
3113c54f92cSHong Zhang      ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3123c54f92cSHong Zhang    }else {
3133c54f92cSHong Zhang      ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3143c54f92cSHong Zhang    }
3152ca6e920SHong Zhang    ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
3162ca6e920SHong Zhang 
3172ca6e920SHong Zhang    /* Solve LHS X = RHS */
318abc2977eSBarry Smith    for (nadj=0; nadj<ts->numcost; nadj++) {
319b5e0532dSHong Zhang      ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
3202ca6e920SHong Zhang    }
3213c54f92cSHong Zhang 
32236eaed60SHong Zhang    /* Update sensitivities, and evaluate integrals if there is any */
323b5e0532dSHong Zhang    if(th->endpoint) { /* two-stage case */
324b5e0532dSHong Zhang      if (th->Theta!=1.) {
3253fd52205SHong Zhang        shift = -1./((th->Theta-1.)*ts->time_step);
326b5e0532dSHong Zhang        ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3272c39e106SBarry Smith        if (ts->vec_costintegral) {
328b5e0532dSHong Zhang          ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
3298263dbbfSHong Zhang        }
330abc2977eSBarry Smith        for (nadj=0; nadj<ts->numcost; nadj++) {
331b5e0532dSHong Zhang          ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
3322c39e106SBarry Smith          if (ts->vec_costintegral) {
33336eaed60SHong Zhang            ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
33436eaed60SHong Zhang          }
3353c54f92cSHong Zhang          ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
3362ca6e920SHong Zhang        }
337b5e0532dSHong Zhang      }else { /* backward Euler */
338b5e0532dSHong Zhang        shift = 0.0;
339b5e0532dSHong Zhang        ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
340b5e0532dSHong Zhang        for (nadj=0; nadj<ts->numcost; nadj++) {
341b5e0532dSHong Zhang          ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
342b5e0532dSHong Zhang          ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
343b5e0532dSHong Zhang          if (ts->vec_costintegral) {
344b5e0532dSHong Zhang            ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
345b5e0532dSHong Zhang          }
346b5e0532dSHong Zhang        }
347b5e0532dSHong Zhang      }
3483fd52205SHong Zhang 
3493fd52205SHong Zhang      if (ts->vecs_sensip) { /* sensitivities wrt parameters */
3505bf8c567SBarry Smith        ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr);
351abc2977eSBarry Smith        for (nadj=0; nadj<ts->numcost; nadj++) {
352b5e0532dSHong Zhang          ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
353b5e0532dSHong Zhang          ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
3543fd52205SHong Zhang        }
355b5e0532dSHong Zhang        if (th->Theta!=1.) {
356b5e0532dSHong Zhang          ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr);
357abc2977eSBarry Smith          for (nadj=0; nadj<ts->numcost; nadj++) {
358b5e0532dSHong Zhang            ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
359b5e0532dSHong Zhang            ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
360b5e0532dSHong Zhang          }
3613fd52205SHong Zhang        }
3622c39e106SBarry Smith        if (ts->vec_costintegral) {
363d4aa0a58SBarry Smith          ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
364abc2977eSBarry Smith          for (nadj=0; nadj<ts->numcost; nadj++) {
36536eaed60SHong Zhang            ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
36636eaed60SHong Zhang          }
367b5e0532dSHong Zhang          if (th->Theta!=1.) {
368b5e0532dSHong Zhang            ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
369abc2977eSBarry Smith            for (nadj=0; nadj<ts->numcost; nadj++) {
37036eaed60SHong Zhang              ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr);
37136eaed60SHong Zhang            }
37236eaed60SHong Zhang          }
3733fd52205SHong Zhang        }
374b5e0532dSHong Zhang      }
3753fd52205SHong Zhang    }else { /* one-stage case */
3763c54f92cSHong Zhang      shift = 0.0;
377a2ae3dd2SHong Zhang      ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
3782c39e106SBarry Smith      if (ts->vec_costintegral) {
379d4aa0a58SBarry Smith        ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
3803c54f92cSHong Zhang      }
381abc2977eSBarry Smith      for (nadj=0; nadj<ts->numcost; nadj++) {
382b5e0532dSHong Zhang        ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
383b5e0532dSHong Zhang        ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
3842c39e106SBarry Smith        if (ts->vec_costintegral) {
38536eaed60SHong Zhang          ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
38636eaed60SHong Zhang        }
3872ca6e920SHong Zhang      }
3883fd52205SHong Zhang      if (ts->vecs_sensip) {
3895bf8c567SBarry Smith        ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
390abc2977eSBarry Smith        for (nadj=0; nadj<ts->numcost; nadj++) {
391b5e0532dSHong Zhang          ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
392b5e0532dSHong Zhang          ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
3933fd52205SHong Zhang        }
3942c39e106SBarry Smith        if (ts->vec_costintegral) {
395d4aa0a58SBarry Smith          ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
396abc2977eSBarry Smith          for (nadj=0; nadj<ts->numcost; nadj++) {
39736eaed60SHong Zhang            ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
39836eaed60SHong Zhang          }
39936eaed60SHong Zhang        }
4003fd52205SHong Zhang      }
4012ca6e920SHong Zhang    }
4022ca6e920SHong Zhang 
4032ca6e920SHong Zhang    th->status = TS_STEP_COMPLETE;
4042ca6e920SHong Zhang    PetscFunctionReturn(0);
4052ca6e920SHong Zhang  }
4062ca6e920SHong Zhang 
4072ca6e920SHong Zhang  #undef __FUNCT__
408cd652676SJed Brown  #define __FUNCT__ "TSInterpolate_Theta"
409cd652676SJed Brown  static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
410cd652676SJed Brown  {
411cd652676SJed Brown    TS_Theta       *th = (TS_Theta*)ts->data;
412be5899b3SLisandro Dalcin    PetscReal      dt  = t - ts->ptime;
413cd652676SJed Brown    PetscErrorCode ierr;
414cd652676SJed Brown 
415cd652676SJed Brown    PetscFunctionBegin;
416a43b19c4SJed Brown    ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
417be5899b3SLisandro Dalcin    if (th->endpoint) dt *= th->Theta;
418be5899b3SLisandro Dalcin    ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
419cd652676SJed Brown    PetscFunctionReturn(0);
420cd652676SJed Brown  }
421cd652676SJed Brown 
4221566a47fSLisandro Dalcin #undef __FUNCT__
4231566a47fSLisandro Dalcin #define __FUNCT__ "TSEvaluateWLTE_Theta"
4241566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
4251566a47fSLisandro Dalcin {
4261566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
4271566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
4281566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
4291566a47fSLisandro Dalcin   PetscErrorCode ierr;
4301566a47fSLisandro Dalcin 
4311566a47fSLisandro Dalcin   PetscFunctionBegin;
4321566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
433*fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
4341566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
435*fecfb714SLisandro Dalcin   {
436be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
437be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
4381566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
4391566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
4401566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
4411566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
4421566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
4431566a47fSLisandro Dalcin     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte);CHKERRQ(ierr);
4441566a47fSLisandro Dalcin   }
4451566a47fSLisandro Dalcin   if (order) *order = 2;
4461566a47fSLisandro Dalcin   PetscFunctionReturn(0);
4471566a47fSLisandro Dalcin }
4481566a47fSLisandro Dalcin 
4491566a47fSLisandro Dalcin #undef __FUNCT__
4501566a47fSLisandro Dalcin #define __FUNCT__ "TSRollBack_Theta"
4511566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
4521566a47fSLisandro Dalcin {
4531566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
4541566a47fSLisandro Dalcin   PetscErrorCode ierr;
4551566a47fSLisandro Dalcin 
4561566a47fSLisandro Dalcin   PetscFunctionBegin;
4571566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
4581566a47fSLisandro Dalcin   if (ts->vec_costintegral && ts->costintegralfwd) {
4591566a47fSLisandro Dalcin     ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
4601566a47fSLisandro Dalcin   }
4611566a47fSLisandro Dalcin   PetscFunctionReturn(0);
4621566a47fSLisandro Dalcin }
4631566a47fSLisandro Dalcin 
464316643e7SJed Brown /*------------------------------------------------------------*/
465316643e7SJed Brown #undef __FUNCT__
466277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
467277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
468316643e7SJed Brown {
469316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
470316643e7SJed Brown   PetscErrorCode ierr;
471316643e7SJed Brown 
472316643e7SJed Brown   PetscFunctionBegin;
4736bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
4746bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
4753b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
476eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
4771566a47fSLisandro Dalcin 
4781566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
4791566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
4801566a47fSLisandro Dalcin 
481b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
482b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
483b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
484b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
485277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
486277b19d0SLisandro Dalcin }
487277b19d0SLisandro Dalcin 
488277b19d0SLisandro Dalcin #undef __FUNCT__
489277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
490277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
491277b19d0SLisandro Dalcin {
492277b19d0SLisandro Dalcin   PetscErrorCode ierr;
493277b19d0SLisandro Dalcin 
494277b19d0SLisandro Dalcin   PetscFunctionBegin;
495277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
496277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
497bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
498bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
499bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
500bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
501316643e7SJed Brown   PetscFunctionReturn(0);
502316643e7SJed Brown }
503316643e7SJed Brown 
504316643e7SJed Brown /*
505316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
5062b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
507316643e7SJed Brown */
508316643e7SJed Brown #undef __FUNCT__
5090f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
5100f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
511316643e7SJed Brown {
512316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
513316643e7SJed Brown   PetscErrorCode ierr;
5147445fe48SJed Brown   Vec            X0,Xdot;
5157445fe48SJed Brown   DM             dm,dmsave;
5161566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
517316643e7SJed Brown 
518316643e7SJed Brown   PetscFunctionBegin;
5197445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
5205a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
5217445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
522b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
5237445fe48SJed Brown 
5247445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
5257445fe48SJed Brown   dmsave = ts->dm;
5267445fe48SJed Brown   ts->dm = dm;
5277445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
5287445fe48SJed Brown   ts->dm = dmsave;
5290d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
530316643e7SJed Brown   PetscFunctionReturn(0);
531316643e7SJed Brown }
532316643e7SJed Brown 
533316643e7SJed Brown #undef __FUNCT__
5340f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
535d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
536316643e7SJed Brown {
537316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
538316643e7SJed Brown   PetscErrorCode ierr;
5397445fe48SJed Brown   Vec            Xdot;
5407445fe48SJed Brown   DM             dm,dmsave;
5411566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
542316643e7SJed Brown 
543316643e7SJed Brown   PetscFunctionBegin;
5447445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
545be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
5460298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
5477445fe48SJed Brown 
5487445fe48SJed Brown   dmsave = ts->dm;
5497445fe48SJed Brown   ts->dm = dm;
550d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
5517445fe48SJed Brown   ts->dm = dmsave;
5520298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
553316643e7SJed Brown   PetscFunctionReturn(0);
554316643e7SJed Brown }
555316643e7SJed Brown 
556316643e7SJed Brown #undef __FUNCT__
557316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
558316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
559316643e7SJed Brown {
560316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
561316643e7SJed Brown   PetscErrorCode ierr;
562316643e7SJed Brown 
563316643e7SJed Brown   PetscFunctionBegin;
564b1cb36f3SHong Zhang   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
565b1cb36f3SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
566b1cb36f3SHong Zhang   }
56739d13666SHong Zhang   if (!th->X) {
568316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
56939d13666SHong Zhang   }
57039d13666SHong Zhang   if (!th->Xdot) {
571316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
57239d13666SHong Zhang   }
57339d13666SHong Zhang   if (!th->X0) {
5743b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
57539d13666SHong Zhang   }
5761566a47fSLisandro Dalcin   if (th->endpoint) {
5771566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
5787445fe48SJed Brown   }
5793b1890cdSShri Abhyankar 
5801566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
5811566a47fSLisandro Dalcin 
5821566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
5831566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
5841566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
5851566a47fSLisandro Dalcin 
5861566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
5871566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
588ef749922SLisandro Dalcin   if (!th->adapt) {
5891566a47fSLisandro Dalcin     ierr = TSAdaptSetType(ts->adapt,TSADAPTNONE);CHKERRQ(ierr);
5901566a47fSLisandro Dalcin   } else {
5911566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
5921566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
5931566a47fSLisandro Dalcin     if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED)
5941566a47fSLisandro Dalcin       ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
5953b1890cdSShri Abhyankar   }
5961566a47fSLisandro Dalcin 
5971566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
598b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
599b4dd3bf9SBarry Smith }
6000c7376e5SHong Zhang 
601b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
602b4dd3bf9SBarry Smith 
603b4dd3bf9SBarry Smith #undef __FUNCT__
604b4dd3bf9SBarry Smith #define __FUNCT__ "TSAdjointSetUp_Theta"
605b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
606b4dd3bf9SBarry Smith {
607b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
608b4dd3bf9SBarry Smith   PetscErrorCode ierr;
609b4dd3bf9SBarry Smith 
610b4dd3bf9SBarry Smith   PetscFunctionBegin;
611b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
6122ca6e920SHong Zhang   if(ts->vecs_sensip) {
613b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
6142ca6e920SHong Zhang   }
615b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
616316643e7SJed Brown   PetscFunctionReturn(0);
617316643e7SJed Brown }
618316643e7SJed Brown /*------------------------------------------------------------*/
619316643e7SJed Brown 
620316643e7SJed Brown #undef __FUNCT__
621316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
6224416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
623316643e7SJed Brown {
624316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
625316643e7SJed Brown   PetscErrorCode ierr;
626316643e7SJed Brown 
627316643e7SJed Brown   PetscFunctionBegin;
628e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
629316643e7SJed Brown   {
6300298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
6310298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
6320298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
63303842d09SLisandro Dalcin     ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
634316643e7SJed Brown   }
635316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
636316643e7SJed Brown   PetscFunctionReturn(0);
637316643e7SJed Brown }
638316643e7SJed Brown 
639316643e7SJed Brown #undef __FUNCT__
640316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
641316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
642316643e7SJed Brown {
643316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
644ace3abfcSBarry Smith   PetscBool      iascii;
645316643e7SJed Brown   PetscErrorCode ierr;
646316643e7SJed Brown 
647316643e7SJed Brown   PetscFunctionBegin;
648251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
649316643e7SJed Brown   if (iascii) {
6507c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
651ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
652316643e7SJed Brown   }
6539c334d8fSLisandro Dalcin   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
654ac75fa18SLisandro Dalcin   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
655316643e7SJed Brown   PetscFunctionReturn(0);
656316643e7SJed Brown }
657316643e7SJed Brown 
6580de4c49aSJed Brown #undef __FUNCT__
6590de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
660560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
6610de4c49aSJed Brown {
6620de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
6630de4c49aSJed Brown 
6640de4c49aSJed Brown   PetscFunctionBegin;
6650de4c49aSJed Brown   *theta = th->Theta;
6660de4c49aSJed Brown   PetscFunctionReturn(0);
6670de4c49aSJed Brown }
6680de4c49aSJed Brown 
6690de4c49aSJed Brown #undef __FUNCT__
6700de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
671560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
6720de4c49aSJed Brown {
6730de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
6740de4c49aSJed Brown 
6750de4c49aSJed Brown   PetscFunctionBegin;
6767c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
6770de4c49aSJed Brown   th->Theta = theta;
6781566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
6790de4c49aSJed Brown   PetscFunctionReturn(0);
6800de4c49aSJed Brown }
681eb284becSJed Brown 
682eb284becSJed Brown #undef __FUNCT__
68378e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
684560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
68526f2ff8fSLisandro Dalcin {
68626f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
68726f2ff8fSLisandro Dalcin 
68826f2ff8fSLisandro Dalcin   PetscFunctionBegin;
68926f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
69026f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
69126f2ff8fSLisandro Dalcin }
69226f2ff8fSLisandro Dalcin 
69326f2ff8fSLisandro Dalcin #undef __FUNCT__
69426f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
695560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
696eb284becSJed Brown {
697eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
698eb284becSJed Brown 
699eb284becSJed Brown   PetscFunctionBegin;
700eb284becSJed Brown   th->endpoint = flg;
701eb284becSJed Brown   PetscFunctionReturn(0);
702eb284becSJed Brown }
7030de4c49aSJed Brown 
704f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
705f9c1d6abSBarry Smith #undef __FUNCT__
706f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
707f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
708f9c1d6abSBarry Smith {
709f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
710f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
7113fd8ae06SJed Brown   const PetscReal one = 1.0;
712f9c1d6abSBarry Smith 
713f9c1d6abSBarry Smith   PetscFunctionBegin;
7143fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
715f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
716f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
717f9c1d6abSBarry Smith   PetscFunctionReturn(0);
718f9c1d6abSBarry Smith }
719f9c1d6abSBarry Smith #endif
720f9c1d6abSBarry Smith 
72142682096SHong Zhang #undef __FUNCT__
72242682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta"
72342682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
72442682096SHong Zhang {
72542682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
72642682096SHong Zhang 
72742682096SHong Zhang   PetscFunctionBegin;
7281566a47fSLisandro Dalcin   if (ns) *ns = 1;
7291566a47fSLisandro Dalcin   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
73042682096SHong Zhang   PetscFunctionReturn(0);
73142682096SHong Zhang }
732f9c1d6abSBarry Smith 
733316643e7SJed Brown /* ------------------------------------------------------------ */
734316643e7SJed Brown /*MC
73596f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
736316643e7SJed Brown 
737316643e7SJed Brown    Level: beginner
738316643e7SJed Brown 
7394eb428fdSBarry Smith    Options Database:
740c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
741c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
74203842d09SLisandro Dalcin .  -ts_theta_adapt <flg> - Use time-step adaptivity with the Theta method
74303842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
7444eb428fdSBarry Smith 
745eb284becSJed Brown    Notes:
7460c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
7470c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
7484eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
7494eb428fdSBarry Smith 
750eb284becSJed Brown    This method can be applied to DAE.
751eb284becSJed Brown 
752eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
753eb284becSJed Brown 
754eb284becSJed Brown .vb
755eb284becSJed Brown   Theta | Theta
756eb284becSJed Brown   -------------
757eb284becSJed Brown         |  1
758eb284becSJed Brown .ve
759eb284becSJed Brown 
760eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
761eb284becSJed Brown 
762eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
763eb284becSJed Brown 
764eb284becSJed Brown .vb
765eb284becSJed Brown   0 | 0         0
766eb284becSJed Brown   1 | 1-Theta   Theta
767eb284becSJed Brown   -------------------
768eb284becSJed Brown     | 1-Theta   Theta
769eb284becSJed Brown .ve
770eb284becSJed Brown 
771eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
772eb284becSJed Brown 
773eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
774eb284becSJed Brown 
775eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
776eb284becSJed Brown 
7774eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
778eb284becSJed Brown 
779eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
780316643e7SJed Brown 
781316643e7SJed Brown M*/
782316643e7SJed Brown #undef __FUNCT__
783316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
7848cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
785316643e7SJed Brown {
786316643e7SJed Brown   TS_Theta       *th;
787316643e7SJed Brown   PetscErrorCode ierr;
788316643e7SJed Brown 
789316643e7SJed Brown   PetscFunctionBegin;
790277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
791316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
792316643e7SJed Brown   ts->ops->view            = TSView_Theta;
793316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
79442f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
795316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
796cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
7971566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
79824655328SShri   ts->ops->rollback        = TSRollBack_Theta;
799316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
8000f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
8010f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
802f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
803f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
804f9c1d6abSBarry Smith #endif
80542682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
80642f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
807b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
808b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
809316643e7SJed Brown 
810b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
811316643e7SJed Brown   ts->data = (void*)th;
812316643e7SJed Brown 
8136f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
814316643e7SJed Brown   th->Theta       = 0.5;
8151566a47fSLisandro Dalcin   th->order       = 2;
8163b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
817bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
818bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
819bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
820bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
821316643e7SJed Brown   PetscFunctionReturn(0);
822316643e7SJed Brown }
8230de4c49aSJed Brown 
8240de4c49aSJed Brown #undef __FUNCT__
8250de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
8260de4c49aSJed Brown /*@
8270de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
8280de4c49aSJed Brown 
8290de4c49aSJed Brown   Not Collective
8300de4c49aSJed Brown 
8310de4c49aSJed Brown   Input Parameter:
8320de4c49aSJed Brown .  ts - timestepping context
8330de4c49aSJed Brown 
8340de4c49aSJed Brown   Output Parameter:
8350de4c49aSJed Brown .  theta - stage abscissa
8360de4c49aSJed Brown 
8370de4c49aSJed Brown   Note:
8380de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
8390de4c49aSJed Brown 
8400de4c49aSJed Brown   Level: Advanced
8410de4c49aSJed Brown 
8420de4c49aSJed Brown .seealso: TSThetaSetTheta()
8430de4c49aSJed Brown @*/
8447087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
8450de4c49aSJed Brown {
8464ac538c5SBarry Smith   PetscErrorCode ierr;
8470de4c49aSJed Brown 
8480de4c49aSJed Brown   PetscFunctionBegin;
849afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8500de4c49aSJed Brown   PetscValidPointer(theta,2);
8514ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
8520de4c49aSJed Brown   PetscFunctionReturn(0);
8530de4c49aSJed Brown }
8540de4c49aSJed Brown 
8550de4c49aSJed Brown #undef __FUNCT__
8560de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
8570de4c49aSJed Brown /*@
8580de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
8590de4c49aSJed Brown 
8600de4c49aSJed Brown   Not Collective
8610de4c49aSJed Brown 
8620de4c49aSJed Brown   Input Parameter:
8630de4c49aSJed Brown +  ts - timestepping context
8640de4c49aSJed Brown -  theta - stage abscissa
8650de4c49aSJed Brown 
8660de4c49aSJed Brown   Options Database:
8670de4c49aSJed Brown .  -ts_theta_theta <theta>
8680de4c49aSJed Brown 
8690de4c49aSJed Brown   Level: Intermediate
8700de4c49aSJed Brown 
8710de4c49aSJed Brown .seealso: TSThetaGetTheta()
8720de4c49aSJed Brown @*/
8737087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
8740de4c49aSJed Brown {
8754ac538c5SBarry Smith   PetscErrorCode ierr;
8760de4c49aSJed Brown 
8770de4c49aSJed Brown   PetscFunctionBegin;
878afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8794ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
8800de4c49aSJed Brown   PetscFunctionReturn(0);
8810de4c49aSJed Brown }
882f33bbcb6SJed Brown 
883eb284becSJed Brown #undef __FUNCT__
88426f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
88526f2ff8fSLisandro Dalcin /*@
88626f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
88726f2ff8fSLisandro Dalcin 
88826f2ff8fSLisandro Dalcin   Not Collective
88926f2ff8fSLisandro Dalcin 
89026f2ff8fSLisandro Dalcin   Input Parameter:
89126f2ff8fSLisandro Dalcin .  ts - timestepping context
89226f2ff8fSLisandro Dalcin 
89326f2ff8fSLisandro Dalcin   Output Parameter:
89426f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
89526f2ff8fSLisandro Dalcin 
89626f2ff8fSLisandro Dalcin   Level: Advanced
89726f2ff8fSLisandro Dalcin 
89826f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
89926f2ff8fSLisandro Dalcin @*/
90026f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
90126f2ff8fSLisandro Dalcin {
90226f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
90326f2ff8fSLisandro Dalcin 
90426f2ff8fSLisandro Dalcin   PetscFunctionBegin;
90526f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
90626f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
907163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
90826f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
90926f2ff8fSLisandro Dalcin }
91026f2ff8fSLisandro Dalcin 
91126f2ff8fSLisandro Dalcin #undef __FUNCT__
912eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
913eb284becSJed Brown /*@
914eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
915eb284becSJed Brown 
916eb284becSJed Brown   Not Collective
917eb284becSJed Brown 
918eb284becSJed Brown   Input Parameter:
919eb284becSJed Brown +  ts - timestepping context
920eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
921eb284becSJed Brown 
922eb284becSJed Brown   Options Database:
923eb284becSJed Brown .  -ts_theta_endpoint <flg>
924eb284becSJed Brown 
925eb284becSJed Brown   Level: Intermediate
926eb284becSJed Brown 
927eb284becSJed Brown .seealso: TSTHETA, TSCN
928eb284becSJed Brown @*/
929eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
930eb284becSJed Brown {
931eb284becSJed Brown   PetscErrorCode ierr;
932eb284becSJed Brown 
933eb284becSJed Brown   PetscFunctionBegin;
934eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
935eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
936eb284becSJed Brown   PetscFunctionReturn(0);
937eb284becSJed Brown }
938eb284becSJed Brown 
939f33bbcb6SJed Brown /*
940f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
941f33bbcb6SJed Brown  * The creation functions for these specializations are below.
942f33bbcb6SJed Brown  */
943f33bbcb6SJed Brown 
944f33bbcb6SJed Brown #undef __FUNCT__
9451566a47fSLisandro Dalcin #define __FUNCT__ "TSSetUp_BEuler"
9461566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
9471566a47fSLisandro Dalcin {
9481566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
9491566a47fSLisandro Dalcin   PetscErrorCode ierr;
9501566a47fSLisandro Dalcin 
9511566a47fSLisandro Dalcin   PetscFunctionBegin;
9521566a47fSLisandro Dalcin   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");
9531566a47fSLisandro Dalcin   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");
9541566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
9551566a47fSLisandro Dalcin   PetscFunctionReturn(0);
9561566a47fSLisandro Dalcin }
9571566a47fSLisandro Dalcin 
9581566a47fSLisandro Dalcin #undef __FUNCT__
959f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
960f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
961f33bbcb6SJed Brown {
962d52bd9f3SBarry Smith   PetscErrorCode ierr;
963d52bd9f3SBarry Smith 
964f33bbcb6SJed Brown   PetscFunctionBegin;
9659c334d8fSLisandro Dalcin   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
9661566a47fSLisandro Dalcin   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
967f33bbcb6SJed Brown   PetscFunctionReturn(0);
968f33bbcb6SJed Brown }
969f33bbcb6SJed Brown 
970f33bbcb6SJed Brown /*MC
971f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
972f33bbcb6SJed Brown 
973f33bbcb6SJed Brown   Level: beginner
974f33bbcb6SJed Brown 
9754eb428fdSBarry Smith   Notes:
976c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
9774eb428fdSBarry Smith 
9781566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
9794eb428fdSBarry Smith 
980f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
981f33bbcb6SJed Brown 
982f33bbcb6SJed Brown M*/
983f33bbcb6SJed Brown #undef __FUNCT__
984f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
9858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
986f33bbcb6SJed Brown {
987f33bbcb6SJed Brown   PetscErrorCode ierr;
988f33bbcb6SJed Brown 
989f33bbcb6SJed Brown   PetscFunctionBegin;
990f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
991f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
9921566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
9930c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
994f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
995f33bbcb6SJed Brown   PetscFunctionReturn(0);
996f33bbcb6SJed Brown }
997f33bbcb6SJed Brown 
998f33bbcb6SJed Brown #undef __FUNCT__
9991566a47fSLisandro Dalcin #define __FUNCT__ "TSSetUp_CN"
10001566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
10011566a47fSLisandro Dalcin {
10021566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
10031566a47fSLisandro Dalcin   PetscErrorCode ierr;
10041566a47fSLisandro Dalcin 
10051566a47fSLisandro Dalcin   PetscFunctionBegin;
10061566a47fSLisandro Dalcin   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");
10071566a47fSLisandro Dalcin   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");
10081566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
10091566a47fSLisandro Dalcin   PetscFunctionReturn(0);
10101566a47fSLisandro Dalcin }
10111566a47fSLisandro Dalcin 
10121566a47fSLisandro Dalcin #undef __FUNCT__
1013f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
1014f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1015f33bbcb6SJed Brown {
1016d52bd9f3SBarry Smith   PetscErrorCode ierr;
1017d52bd9f3SBarry Smith 
1018f33bbcb6SJed Brown   PetscFunctionBegin;
10199c334d8fSLisandro Dalcin   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
10201566a47fSLisandro Dalcin   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
1021f33bbcb6SJed Brown   PetscFunctionReturn(0);
1022f33bbcb6SJed Brown }
1023f33bbcb6SJed Brown 
1024f33bbcb6SJed Brown /*MC
1025f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1026f33bbcb6SJed Brown 
1027f33bbcb6SJed Brown   Level: beginner
1028f33bbcb6SJed Brown 
1029f33bbcb6SJed Brown   Notes:
10307cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
10317cf5af47SJed Brown 
10327cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1033f33bbcb6SJed Brown 
1034f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1035f33bbcb6SJed Brown 
1036f33bbcb6SJed Brown M*/
1037f33bbcb6SJed Brown #undef __FUNCT__
1038f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
10398cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1040f33bbcb6SJed Brown {
1041f33bbcb6SJed Brown   PetscErrorCode ierr;
1042f33bbcb6SJed Brown 
1043f33bbcb6SJed Brown   PetscFunctionBegin;
1044f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1045f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1046eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
10470c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1048f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1049f33bbcb6SJed Brown   PetscFunctionReturn(0);
1050f33bbcb6SJed Brown }
1051