xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision c82ed9623bc85b782d7be52242d56e82caddbc18)
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 {
10316643e7SJed Brown   Vec          X,Xdot;                   /* Storage for one stage */
113b1890cdSShri Abhyankar   Vec          X0;                       /* work vector to store X0 */
12eb284becSJed Brown   Vec          affine;                   /* Affine vector needed for residual at beginning of step */
13b5e0532dSHong Zhang   Vec          *VecsDeltaLam;             /* Increment of the adjoint sensitivity w.r.t IC at stage*/
14b5e0532dSHong Zhang   Vec          *VecsDeltaMu;              /* Increment of the adjoint sensitivity w.r.t P at stage*/
15b5e0532dSHong Zhang   Vec          *VecsSensiTemp;            /* Vector to be timed with Jacobian transpose*/
16ace3abfcSBarry Smith   PetscBool    extrapolate;
17eb284becSJed Brown   PetscBool    endpoint;
18316643e7SJed Brown   PetscReal    Theta;
19316643e7SJed Brown   PetscReal    stage_time;
203b1890cdSShri Abhyankar   TSStepStatus status;
213b1890cdSShri Abhyankar   char         *name;
223b1890cdSShri Abhyankar   PetscInt     order;
233b1890cdSShri Abhyankar   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
243b1890cdSShri Abhyankar   PetscBool    adapt;  /* use time-step adaptivity ? */
25b5e0532dSHong Zhang   PetscReal    ptime;
26316643e7SJed Brown } TS_Theta;
27316643e7SJed Brown 
28316643e7SJed Brown #undef __FUNCT__
297445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot"
307445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
317445fe48SJed Brown {
327445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
337445fe48SJed Brown   PetscErrorCode ierr;
347445fe48SJed Brown 
357445fe48SJed Brown   PetscFunctionBegin;
367445fe48SJed Brown   if (X0) {
377445fe48SJed Brown     if (dm && dm != ts->dm) {
380d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
397445fe48SJed Brown     } else *X0 = ts->vec_sol;
407445fe48SJed Brown   }
417445fe48SJed Brown   if (Xdot) {
427445fe48SJed Brown     if (dm && dm != ts->dm) {
430d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
447445fe48SJed Brown     } else *Xdot = th->Xdot;
457445fe48SJed Brown   }
467445fe48SJed Brown   PetscFunctionReturn(0);
477445fe48SJed Brown }
487445fe48SJed Brown 
490d0b770aSPeter Brune #undef __FUNCT__
500d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot"
510d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
520d0b770aSPeter Brune {
530d0b770aSPeter Brune   PetscErrorCode ierr;
540d0b770aSPeter Brune 
550d0b770aSPeter Brune   PetscFunctionBegin;
560d0b770aSPeter Brune   if (X0) {
570d0b770aSPeter Brune     if (dm && dm != ts->dm) {
580d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
590d0b770aSPeter Brune     }
600d0b770aSPeter Brune   }
610d0b770aSPeter Brune   if (Xdot) {
620d0b770aSPeter Brune     if (dm && dm != ts->dm) {
630d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
640d0b770aSPeter Brune     }
650d0b770aSPeter Brune   }
660d0b770aSPeter Brune   PetscFunctionReturn(0);
670d0b770aSPeter Brune }
680d0b770aSPeter Brune 
697445fe48SJed Brown #undef __FUNCT__
707445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta"
717445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
727445fe48SJed Brown {
737445fe48SJed Brown 
747445fe48SJed Brown   PetscFunctionBegin;
757445fe48SJed Brown   PetscFunctionReturn(0);
767445fe48SJed Brown }
777445fe48SJed Brown 
787445fe48SJed Brown #undef __FUNCT__
797445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta"
807445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
817445fe48SJed Brown {
827445fe48SJed Brown   TS             ts = (TS)ctx;
837445fe48SJed Brown   PetscErrorCode ierr;
847445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
857445fe48SJed Brown 
867445fe48SJed Brown   PetscFunctionBegin;
877445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
887445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
897445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
907445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
917445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
927445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
930d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
940d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
957445fe48SJed Brown   PetscFunctionReturn(0);
967445fe48SJed Brown }
977445fe48SJed Brown 
987445fe48SJed Brown #undef __FUNCT__
99258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta"
100258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
101258e1594SPeter Brune {
102258e1594SPeter Brune 
103258e1594SPeter Brune   PetscFunctionBegin;
104258e1594SPeter Brune   PetscFunctionReturn(0);
105258e1594SPeter Brune }
106258e1594SPeter Brune 
107258e1594SPeter Brune #undef __FUNCT__
108258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta"
109258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
110258e1594SPeter Brune {
111258e1594SPeter Brune   TS             ts = (TS)ctx;
112258e1594SPeter Brune   PetscErrorCode ierr;
113258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
114258e1594SPeter Brune 
115258e1594SPeter Brune   PetscFunctionBegin;
116258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
117258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
118258e1594SPeter Brune 
119258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121258e1594SPeter Brune 
122258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124258e1594SPeter Brune 
125258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
126258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
127258e1594SPeter Brune   PetscFunctionReturn(0);
128258e1594SPeter Brune }
129258e1594SPeter Brune 
1303b1890cdSShri Abhyankar #undef __FUNCT__
1313b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta"
1323b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
1333b1890cdSShri Abhyankar {
1343b1890cdSShri Abhyankar   PetscErrorCode ierr;
1353b1890cdSShri Abhyankar   TS_Theta       *th = (TS_Theta*)ts->data;
1363b1890cdSShri Abhyankar 
1373b1890cdSShri Abhyankar   PetscFunctionBegin;
138ce94432eSBarry Smith   if (order == 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"No time-step adaptivity implemented for 1st order theta method; Run with -ts_adapt_type none");
1393b1890cdSShri Abhyankar   if (order == th->order) {
1403b1890cdSShri Abhyankar     if (th->endpoint) {
1413b1890cdSShri Abhyankar       ierr = VecCopy(th->X,U);CHKERRQ(ierr);
1423b1890cdSShri Abhyankar     } else {
1433b1890cdSShri Abhyankar       PetscReal shift = 1./(th->Theta*ts->time_step);
1443b1890cdSShri Abhyankar       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr);
1453b1890cdSShri Abhyankar       ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr);
1463b1890cdSShri Abhyankar     }
1473b1890cdSShri Abhyankar   } else if (order == th->order-1 && order) {
1483b1890cdSShri Abhyankar     ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr);
1493b1890cdSShri Abhyankar   }
1503b1890cdSShri Abhyankar   PetscFunctionReturn(0);
1513b1890cdSShri Abhyankar }
152258e1594SPeter Brune 
153258e1594SPeter Brune #undef __FUNCT__
15424655328SShri #define __FUNCT__ "TSRollBack_Theta"
15524655328SShri static PetscErrorCode TSRollBack_Theta(TS ts)
15624655328SShri {
15724655328SShri   TS_Theta       *th = (TS_Theta*)ts->data;
15824655328SShri   PetscErrorCode ierr;
15924655328SShri 
16024655328SShri   PetscFunctionBegin;
16124655328SShri   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
16224655328SShri   th->status    = TS_STEP_INCOMPLETE;
16324655328SShri   PetscFunctionReturn(0);
16424655328SShri }
16524655328SShri 
16624655328SShri #undef __FUNCT__
167316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
168193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
169316643e7SJed Brown {
170316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1713b1890cdSShri Abhyankar   PetscInt       its,lits,reject,next_scheme;
1723b1890cdSShri Abhyankar   PetscReal      next_time_step;
1733b1890cdSShri Abhyankar   TSAdapt        adapt;
1744957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
175051f2191SLisandro Dalcin   PetscErrorCode ierr;
176316643e7SJed Brown 
177316643e7SJed Brown   PetscFunctionBegin;
1783b1890cdSShri Abhyankar   th->status = TS_STEP_INCOMPLETE;
1793b1890cdSShri Abhyankar   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
180051f2191SLisandro Dalcin   for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
181b296d7d5SJed Brown     PetscReal shift = 1./(th->Theta*ts->time_step);
182eb284becSJed Brown     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
183b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
184b8123daeSJed Brown     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
185316643e7SJed Brown 
186eb284becSJed Brown     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
187eb284becSJed Brown       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
188eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
189eb284becSJed Brown       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
190eb284becSJed Brown       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
191eb284becSJed Brown     }
192ace68cafSJed Brown     if (th->extrapolate) {
193b296d7d5SJed Brown       ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
194ace68cafSJed Brown     } else {
1952b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
196ace68cafSJed Brown     }
197eb284becSJed Brown     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
198316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
199316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
2005ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
201051f2191SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr);
202552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
203b295832fSPierre Barbier de Reuille     ierr = TSAdaptCheckStage(adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
204051f2191SLisandro Dalcin     if (!stageok) {accept = PETSC_FALSE; goto reject_step;}
205051f2191SLisandro Dalcin 
2060298fd71SBarry Smith     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
207051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2083b1890cdSShri Abhyankar     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
209552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2103b1890cdSShri Abhyankar     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
2110298fd71SBarry Smith     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
2123b1890cdSShri Abhyankar     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
213051f2191SLisandro Dalcin     if (!accept) {           /* Roll back the current step */
214051f2191SLisandro Dalcin       ts->ptime += next_time_step; /* This will be undone in rollback */
215051f2191SLisandro Dalcin       th->status = TS_STEP_INCOMPLETE;
216051f2191SLisandro Dalcin       ierr = TSRollBack(ts);CHKERRQ(ierr);
217051f2191SLisandro Dalcin       goto reject_step;
218051f2191SLisandro Dalcin     }
2193b1890cdSShri Abhyankar 
22017145e75SHong Zhang     if (ts->vec_costintegral) {
22117145e75SHong Zhang       /* Evolve ts->vec_costintegral to compute integrals */
22217145e75SHong Zhang       if (th->endpoint) {
22317145e75SHong Zhang         ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
22417145e75SHong Zhang         ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(1.-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
22517145e75SHong Zhang       }
22617145e75SHong Zhang       ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
22717145e75SHong Zhang       if (th->endpoint) {
22817145e75SHong Zhang         ierr = VecAXPY(ts->vec_costintegral,ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
22917145e75SHong Zhang       }else {
23017145e75SHong Zhang         ierr = VecAXPY(ts->vec_costintegral,ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
23117145e75SHong Zhang       }
23217145e75SHong Zhang     }
23317145e75SHong Zhang 
2343b1890cdSShri Abhyankar     /* ignore next_scheme for now */
2352b5a38e1SLisandro Dalcin     ts->ptime    += ts->time_step;
236cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
237316643e7SJed Brown     ts->steps++;
2383b1890cdSShri Abhyankar     th->status = TS_STEP_COMPLETE;
239051f2191SLisandro Dalcin     break;
240051f2191SLisandro Dalcin 
241051f2191SLisandro Dalcin reject_step:
242051f2191SLisandro Dalcin     if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
243051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
244051f2191SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
2453b1890cdSShri Abhyankar     }
246051f2191SLisandro Dalcin     continue;
2473b1890cdSShri Abhyankar   }
248316643e7SJed Brown   PetscFunctionReturn(0);
249316643e7SJed Brown }
250316643e7SJed Brown 
251cd652676SJed Brown #undef __FUNCT__
25242f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_Theta"
25342f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
2542ca6e920SHong Zhang {
2552ca6e920SHong Zhang   TS_Theta            *th = (TS_Theta*)ts->data;
256b5e0532dSHong Zhang   Vec                 *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
2572ca6e920SHong Zhang   PetscInt            nadj;
2582ca6e920SHong Zhang   PetscErrorCode      ierr;
2592ca6e920SHong Zhang   Mat                 J,Jp;
2602ca6e920SHong Zhang   KSP                 ksp;
2612ca6e920SHong Zhang   PetscReal           shift;
2622ca6e920SHong Zhang 
2632ca6e920SHong Zhang   PetscFunctionBegin;
2642ca6e920SHong Zhang 
2652ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
266302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
2672ca6e920SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
268b5e0532dSHong Zhang 
269b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
2703fd52205SHong Zhang   th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
271b5e0532dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
2722ca6e920SHong Zhang 
2732ca6e920SHong Zhang   ierr = TSPreStep(ts);CHKERRQ(ierr);
2742ca6e920SHong Zhang 
275a4cab896SHong Zhang   /* Build RHS */
2762c39e106SBarry Smith   if (ts->vec_costintegral) { /* Cost function has an integral term */
27705755b9cSHong Zhang     if (th->endpoint) {
278d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr);
27905755b9cSHong Zhang     }else {
280d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
28105755b9cSHong Zhang     }
28205755b9cSHong Zhang   }
283abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
284b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
285b5e0532dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
2862c39e106SBarry Smith     if (ts->vec_costintegral) {
287b5e0532dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
28836eaed60SHong Zhang     }
2892ca6e920SHong Zhang   }
2903c54f92cSHong Zhang 
2912ca6e920SHong Zhang   /* Build LHS */
2922ca6e920SHong Zhang   shift = -1./(th->Theta*ts->time_step);
2933c54f92cSHong Zhang   if (th->endpoint) {
2943c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2953c54f92cSHong Zhang   }else {
2963c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2973c54f92cSHong Zhang   }
2982ca6e920SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
2992ca6e920SHong Zhang 
3002ca6e920SHong Zhang   /* Solve LHS X = RHS */
301abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
302b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
3032ca6e920SHong Zhang   }
3043c54f92cSHong Zhang 
30536eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
306b5e0532dSHong Zhang   if(th->endpoint) { /* two-stage case */
307b5e0532dSHong Zhang     if (th->Theta!=1.) {
3083fd52205SHong Zhang       shift = -1./((th->Theta-1.)*ts->time_step);
309b5e0532dSHong Zhang       ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3102c39e106SBarry Smith       if (ts->vec_costintegral) {
311b5e0532dSHong Zhang         ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
312fb07d950SHong Zhang         if (!ts->costintegralfwd) {
3132c39e106SBarry Smith           /* Evolve ts->vec_costintegral to compute integrals */
314d4aa0a58SBarry Smith           ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
3152c39e106SBarry Smith           ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
316b5e0532dSHong Zhang           ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
3172c39e106SBarry Smith           ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1.),ts->vec_costintegrand);CHKERRQ(ierr);
31836eaed60SHong Zhang         }
3198263dbbfSHong Zhang       }
320abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
321b5e0532dSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
3222c39e106SBarry Smith         if (ts->vec_costintegral) {
32336eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
32436eaed60SHong Zhang         }
3253c54f92cSHong Zhang         ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
3262ca6e920SHong Zhang       }
327b5e0532dSHong Zhang     }else { /* backward Euler */
328b5e0532dSHong Zhang       shift = 0.0;
329b5e0532dSHong Zhang       ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
330b5e0532dSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
331b5e0532dSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
332b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
333b5e0532dSHong Zhang         if (ts->vec_costintegral) {
334b5e0532dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
335b5e0532dSHong Zhang           if (!ts->costintegralfwd) {
336b5e0532dSHong Zhang           /* Evolve ts->vec_costintegral to compute integrals */
337b5e0532dSHong Zhang             ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
338b5e0532dSHong Zhang             ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
339b5e0532dSHong Zhang           }
340b5e0532dSHong Zhang         }
341b5e0532dSHong Zhang       }
342b5e0532dSHong Zhang     }
3433fd52205SHong Zhang 
3443fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
3455bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr);
346abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
347b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
348b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
3493fd52205SHong Zhang       }
350b5e0532dSHong Zhang       if (th->Theta!=1.) {
351b5e0532dSHong Zhang         ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr);
352abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
353b5e0532dSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
354b5e0532dSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
355b5e0532dSHong Zhang         }
3563fd52205SHong Zhang       }
3572c39e106SBarry Smith       if (ts->vec_costintegral) {
358d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
359abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
36036eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
36136eaed60SHong Zhang         }
362b5e0532dSHong Zhang         if (th->Theta!=1.) {
363b5e0532dSHong Zhang           ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,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*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr);
36636eaed60SHong Zhang           }
36736eaed60SHong Zhang         }
3683fd52205SHong Zhang       }
369b5e0532dSHong Zhang     }
3703fd52205SHong Zhang   }else { /* one-stage case */
3713c54f92cSHong Zhang     shift = 0.0;
372a2ae3dd2SHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
3732c39e106SBarry Smith     if (ts->vec_costintegral) {
374d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
375fb07d950SHong Zhang       if (!ts->costintegralfwd) {
3762c39e106SBarry Smith         /* Evolve ts->vec_costintegral to compute integrals */
377d4aa0a58SBarry Smith         ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
3782c39e106SBarry Smith         ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
37936eaed60SHong Zhang       }
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   ts->ptime += ts->time_step;
4042ca6e920SHong Zhang   ts->steps++;
4052ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
4062ca6e920SHong Zhang   PetscFunctionReturn(0);
4072ca6e920SHong Zhang }
4082ca6e920SHong Zhang 
4092ca6e920SHong Zhang #undef __FUNCT__
410cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
411cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
412cd652676SJed Brown {
413cd652676SJed Brown   TS_Theta       *th   = (TS_Theta*)ts->data;
4145a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
415cd652676SJed Brown   PetscErrorCode ierr;
416cd652676SJed Brown 
417cd652676SJed Brown   PetscFunctionBegin;
418a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
4195a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
4205a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
421cd652676SJed Brown   PetscFunctionReturn(0);
422cd652676SJed Brown }
423cd652676SJed Brown 
424316643e7SJed Brown /*------------------------------------------------------------*/
425316643e7SJed Brown #undef __FUNCT__
426277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
427277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
428316643e7SJed Brown {
429316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
430316643e7SJed Brown   PetscErrorCode ierr;
431316643e7SJed Brown 
432316643e7SJed Brown   PetscFunctionBegin;
4336bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
4346bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
4353b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
436eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
437b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
438b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
439b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
440277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
441277b19d0SLisandro Dalcin }
442277b19d0SLisandro Dalcin 
443277b19d0SLisandro Dalcin #undef __FUNCT__
444277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
445277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
446277b19d0SLisandro Dalcin {
447277b19d0SLisandro Dalcin   PetscErrorCode ierr;
448277b19d0SLisandro Dalcin 
449277b19d0SLisandro Dalcin   PetscFunctionBegin;
450277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
451277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
452bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
453bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
454bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
455bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
456316643e7SJed Brown   PetscFunctionReturn(0);
457316643e7SJed Brown }
458316643e7SJed Brown 
459316643e7SJed Brown /*
460316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
4612b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
462316643e7SJed Brown */
463316643e7SJed Brown #undef __FUNCT__
4640f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
4650f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
466316643e7SJed Brown {
467316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
468316643e7SJed Brown   PetscErrorCode ierr;
4697445fe48SJed Brown   Vec            X0,Xdot;
4707445fe48SJed Brown   DM             dm,dmsave;
471b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
472316643e7SJed Brown 
473316643e7SJed Brown   PetscFunctionBegin;
4747445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
4755a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
4767445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
477b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
4787445fe48SJed Brown 
4797445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
4807445fe48SJed Brown   dmsave = ts->dm;
4817445fe48SJed Brown   ts->dm = dm;
4827445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
4837445fe48SJed Brown   ts->dm = dmsave;
4840d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
485316643e7SJed Brown   PetscFunctionReturn(0);
486316643e7SJed Brown }
487316643e7SJed Brown 
488316643e7SJed Brown #undef __FUNCT__
4890f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
490d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
491316643e7SJed Brown {
492316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
493316643e7SJed Brown   PetscErrorCode ierr;
4947445fe48SJed Brown   Vec            Xdot;
4957445fe48SJed Brown   DM             dm,dmsave;
496b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
497316643e7SJed Brown 
498316643e7SJed Brown   PetscFunctionBegin;
4997445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
5007445fe48SJed Brown 
5010f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
5020298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
5037445fe48SJed Brown 
5047445fe48SJed Brown   dmsave = ts->dm;
5057445fe48SJed Brown   ts->dm = dm;
506d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
5077445fe48SJed Brown   ts->dm = dmsave;
5080298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
509316643e7SJed Brown   PetscFunctionReturn(0);
510316643e7SJed Brown }
511316643e7SJed Brown 
512316643e7SJed Brown #undef __FUNCT__
513316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
514316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
515316643e7SJed Brown {
516316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
517316643e7SJed Brown   PetscErrorCode ierr;
5187445fe48SJed Brown   SNES           snes;
519ef749922SLisandro Dalcin   TSAdapt        adapt;
5207445fe48SJed Brown   DM             dm;
521316643e7SJed Brown 
522316643e7SJed Brown   PetscFunctionBegin;
52339d13666SHong Zhang   if (!th->X) {
524316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
52539d13666SHong Zhang   }
52639d13666SHong Zhang   if (!th->Xdot) {
527316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
52839d13666SHong Zhang   }
52939d13666SHong Zhang   if (!th->X0) {
5303b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
53139d13666SHong Zhang   }
5327445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
5337445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
5347445fe48SJed Brown   if (dm) {
5357445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
536258e1594SPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
5377445fe48SJed Brown   }
5383b1890cdSShri Abhyankar   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
5393b1890cdSShri Abhyankar   else th->order = 1;
5403b1890cdSShri Abhyankar 
541552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
542ef749922SLisandro Dalcin   if (!th->adapt) {
5433b1890cdSShri Abhyankar     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
5443b1890cdSShri Abhyankar   }
545b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
546b4dd3bf9SBarry Smith }
5470c7376e5SHong Zhang 
5480c7376e5SHong Zhang #undef __FUNCT__
5490c7376e5SHong Zhang #define __FUNCT__ "TSSetUp_BEuler"
5500c7376e5SHong Zhang static PetscErrorCode TSSetUp_BEuler(TS ts)
5510c7376e5SHong Zhang {
5520c7376e5SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
5530c7376e5SHong Zhang   PetscErrorCode ierr;
5540c7376e5SHong Zhang 
5550c7376e5SHong Zhang   PetscFunctionBegin;
556b9eb5ee8SHong Zhang   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");
557421bc905SSatish Balay   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
5580c7376e5SHong Zhang   PetscFunctionReturn(0);
5590c7376e5SHong Zhang }
5600c7376e5SHong Zhang 
5610c7376e5SHong Zhang #undef __FUNCT__
5620c7376e5SHong Zhang #define __FUNCT__ "TSSetUp_CN"
5630c7376e5SHong Zhang static PetscErrorCode TSSetUp_CN(TS ts)
5640c7376e5SHong Zhang {
5650c7376e5SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
5660c7376e5SHong Zhang   PetscErrorCode ierr;
5670c7376e5SHong Zhang 
5680c7376e5SHong Zhang   PetscFunctionBegin;
569b9eb5ee8SHong Zhang   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");
570b9eb5ee8SHong Zhang   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");
571421bc905SSatish Balay   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
5720c7376e5SHong Zhang   PetscFunctionReturn(0);
5730c7376e5SHong Zhang }
574b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
575b4dd3bf9SBarry Smith 
576b4dd3bf9SBarry Smith #undef __FUNCT__
577b4dd3bf9SBarry Smith #define __FUNCT__ "TSAdjointSetUp_Theta"
578b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
579b4dd3bf9SBarry Smith {
580b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
581b4dd3bf9SBarry Smith   PetscErrorCode ierr;
582b4dd3bf9SBarry Smith 
583b4dd3bf9SBarry Smith   PetscFunctionBegin;
584b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
5852ca6e920SHong Zhang   if(ts->vecs_sensip) {
586b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
5872ca6e920SHong Zhang   }
588b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
589316643e7SJed Brown   PetscFunctionReturn(0);
590316643e7SJed Brown }
591316643e7SJed Brown /*------------------------------------------------------------*/
592316643e7SJed Brown 
593316643e7SJed Brown #undef __FUNCT__
594316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
5958c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts)
596316643e7SJed Brown {
597316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
598316643e7SJed Brown   PetscErrorCode ierr;
599316643e7SJed Brown 
600316643e7SJed Brown   PetscFunctionBegin;
601e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
602316643e7SJed Brown   {
6030298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
6040298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
6050298fd71SBarry 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);
6060298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
607316643e7SJed Brown   }
608316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
609316643e7SJed Brown   PetscFunctionReturn(0);
610316643e7SJed Brown }
611316643e7SJed Brown 
612316643e7SJed Brown #undef __FUNCT__
613316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
614316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
615316643e7SJed Brown {
616316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
617ace3abfcSBarry Smith   PetscBool      iascii;
618316643e7SJed Brown   PetscErrorCode ierr;
619316643e7SJed Brown 
620316643e7SJed Brown   PetscFunctionBegin;
621251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
622316643e7SJed Brown   if (iascii) {
6237c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
624ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
625316643e7SJed Brown   }
626ac75fa18SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
627316643e7SJed Brown   PetscFunctionReturn(0);
628316643e7SJed Brown }
629316643e7SJed Brown 
6300de4c49aSJed Brown #undef __FUNCT__
6310de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
6327087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
6330de4c49aSJed Brown {
6340de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
6350de4c49aSJed Brown 
6360de4c49aSJed Brown   PetscFunctionBegin;
6370de4c49aSJed Brown   *theta = th->Theta;
6380de4c49aSJed Brown   PetscFunctionReturn(0);
6390de4c49aSJed Brown }
6400de4c49aSJed Brown 
6410de4c49aSJed Brown #undef __FUNCT__
6420de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
6437087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
6440de4c49aSJed Brown {
6450de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
6460de4c49aSJed Brown 
6470de4c49aSJed Brown   PetscFunctionBegin;
6487c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
6490de4c49aSJed Brown   th->Theta = theta;
6500de4c49aSJed Brown   PetscFunctionReturn(0);
6510de4c49aSJed Brown }
652eb284becSJed Brown 
653eb284becSJed Brown #undef __FUNCT__
65478e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
65526f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
65626f2ff8fSLisandro Dalcin {
65726f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
65826f2ff8fSLisandro Dalcin 
65926f2ff8fSLisandro Dalcin   PetscFunctionBegin;
66026f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
66126f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
66226f2ff8fSLisandro Dalcin }
66326f2ff8fSLisandro Dalcin 
66426f2ff8fSLisandro Dalcin #undef __FUNCT__
66526f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
666eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
667eb284becSJed Brown {
668eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
669eb284becSJed Brown 
670eb284becSJed Brown   PetscFunctionBegin;
671eb284becSJed Brown   th->endpoint = flg;
672eb284becSJed Brown   PetscFunctionReturn(0);
673eb284becSJed Brown }
6740de4c49aSJed Brown 
675f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
676f9c1d6abSBarry Smith #undef __FUNCT__
677f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
678f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
679f9c1d6abSBarry Smith {
680f9c1d6abSBarry Smith   PetscComplex z   = xr + xi*PETSC_i,f;
681f9c1d6abSBarry Smith   TS_Theta     *th = (TS_Theta*)ts->data;
6823fd8ae06SJed Brown   const PetscReal one = 1.0;
683f9c1d6abSBarry Smith 
684f9c1d6abSBarry Smith   PetscFunctionBegin;
6853fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
686f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
687f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
688f9c1d6abSBarry Smith   PetscFunctionReturn(0);
689f9c1d6abSBarry Smith }
690f9c1d6abSBarry Smith #endif
691f9c1d6abSBarry Smith 
69242682096SHong Zhang #undef __FUNCT__
69342682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta"
69442682096SHong Zhang static PetscErrorCode  TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
69542682096SHong Zhang {
69642682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
69742682096SHong Zhang 
69842682096SHong Zhang   PetscFunctionBegin;
6992ca6e920SHong Zhang   *ns = 1;
7002ca6e920SHong Zhang   if(Y) {
701b5e0532dSHong Zhang     *Y  = (th->endpoint)?&(th->X0):&(th->X);
7022ca6e920SHong Zhang   }
70342682096SHong Zhang   PetscFunctionReturn(0);
70442682096SHong Zhang }
705f9c1d6abSBarry Smith 
706316643e7SJed Brown /* ------------------------------------------------------------ */
707316643e7SJed Brown /*MC
70896f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
709316643e7SJed Brown 
710316643e7SJed Brown    Level: beginner
711316643e7SJed Brown 
7124eb428fdSBarry Smith    Options Database:
713*c82ed962SBarry Smith +      -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
714*c82ed962SBarry Smith .      -ts_theta_extrapolate <flg> - Extrapolate stage solution from previous solution (sometimes unstable)
715*c82ed962SBarry Smith .      -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
716*c82ed962SBarry Smith -     -ts_theta_adapt <flg> - Use time-step adaptivity with the Theta method
7174eb428fdSBarry Smith 
718eb284becSJed Brown    Notes:
7190c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
7200c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
7214eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
7224eb428fdSBarry Smith 
723eb284becSJed Brown    This method can be applied to DAE.
724eb284becSJed Brown 
725eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
726eb284becSJed Brown 
727eb284becSJed Brown .vb
728eb284becSJed Brown   Theta | Theta
729eb284becSJed Brown   -------------
730eb284becSJed Brown         |  1
731eb284becSJed Brown .ve
732eb284becSJed Brown 
733eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
734eb284becSJed Brown 
735eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
736eb284becSJed Brown 
737eb284becSJed Brown .vb
738eb284becSJed Brown   0 | 0         0
739eb284becSJed Brown   1 | 1-Theta   Theta
740eb284becSJed Brown   -------------------
741eb284becSJed Brown     | 1-Theta   Theta
742eb284becSJed Brown .ve
743eb284becSJed Brown 
744eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
745eb284becSJed Brown 
746eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
747eb284becSJed Brown 
748eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
749eb284becSJed Brown 
7504eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
751eb284becSJed Brown 
752eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
753316643e7SJed Brown 
754316643e7SJed Brown M*/
755316643e7SJed Brown #undef __FUNCT__
756316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
7578cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
758316643e7SJed Brown {
759316643e7SJed Brown   TS_Theta       *th;
760316643e7SJed Brown   PetscErrorCode ierr;
761316643e7SJed Brown 
762316643e7SJed Brown   PetscFunctionBegin;
763277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
764316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
765316643e7SJed Brown   ts->ops->view            = TSView_Theta;
766316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
76742f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
768316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
769cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
7703b1890cdSShri Abhyankar   ts->ops->evaluatestep    = TSEvaluateStep_Theta;
77124655328SShri   ts->ops->rollback        = TSRollBack_Theta;
772316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
7730f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
7740f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
775f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
776f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
777f9c1d6abSBarry Smith #endif
77842682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
77942f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
780316643e7SJed Brown 
781b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
782316643e7SJed Brown   ts->data = (void*)th;
783316643e7SJed Brown 
7846f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
785316643e7SJed Brown   th->Theta       = 0.5;
7863b1890cdSShri Abhyankar   th->ccfl        = 1.0;
7873b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
788bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
789bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
790bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
791bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
792316643e7SJed Brown   PetscFunctionReturn(0);
793316643e7SJed Brown }
7940de4c49aSJed Brown 
7950de4c49aSJed Brown #undef __FUNCT__
7960de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
7970de4c49aSJed Brown /*@
7980de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
7990de4c49aSJed Brown 
8000de4c49aSJed Brown   Not Collective
8010de4c49aSJed Brown 
8020de4c49aSJed Brown   Input Parameter:
8030de4c49aSJed Brown .  ts - timestepping context
8040de4c49aSJed Brown 
8050de4c49aSJed Brown   Output Parameter:
8060de4c49aSJed Brown .  theta - stage abscissa
8070de4c49aSJed Brown 
8080de4c49aSJed Brown   Note:
8090de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
8100de4c49aSJed Brown 
8110de4c49aSJed Brown   Level: Advanced
8120de4c49aSJed Brown 
8130de4c49aSJed Brown .seealso: TSThetaSetTheta()
8140de4c49aSJed Brown @*/
8157087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
8160de4c49aSJed Brown {
8174ac538c5SBarry Smith   PetscErrorCode ierr;
8180de4c49aSJed Brown 
8190de4c49aSJed Brown   PetscFunctionBegin;
820afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8210de4c49aSJed Brown   PetscValidPointer(theta,2);
8224ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
8230de4c49aSJed Brown   PetscFunctionReturn(0);
8240de4c49aSJed Brown }
8250de4c49aSJed Brown 
8260de4c49aSJed Brown #undef __FUNCT__
8270de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
8280de4c49aSJed Brown /*@
8290de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
8300de4c49aSJed Brown 
8310de4c49aSJed Brown   Not Collective
8320de4c49aSJed Brown 
8330de4c49aSJed Brown   Input Parameter:
8340de4c49aSJed Brown +  ts - timestepping context
8350de4c49aSJed Brown -  theta - stage abscissa
8360de4c49aSJed Brown 
8370de4c49aSJed Brown   Options Database:
8380de4c49aSJed Brown .  -ts_theta_theta <theta>
8390de4c49aSJed Brown 
8400de4c49aSJed Brown   Level: Intermediate
8410de4c49aSJed Brown 
8420de4c49aSJed Brown .seealso: TSThetaGetTheta()
8430de4c49aSJed Brown @*/
8447087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
8450de4c49aSJed Brown {
8464ac538c5SBarry Smith   PetscErrorCode ierr;
8470de4c49aSJed Brown 
8480de4c49aSJed Brown   PetscFunctionBegin;
849afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8504ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
8510de4c49aSJed Brown   PetscFunctionReturn(0);
8520de4c49aSJed Brown }
853f33bbcb6SJed Brown 
854eb284becSJed Brown #undef __FUNCT__
85526f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
85626f2ff8fSLisandro Dalcin /*@
85726f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
85826f2ff8fSLisandro Dalcin 
85926f2ff8fSLisandro Dalcin   Not Collective
86026f2ff8fSLisandro Dalcin 
86126f2ff8fSLisandro Dalcin   Input Parameter:
86226f2ff8fSLisandro Dalcin .  ts - timestepping context
86326f2ff8fSLisandro Dalcin 
86426f2ff8fSLisandro Dalcin   Output Parameter:
86526f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
86626f2ff8fSLisandro Dalcin 
86726f2ff8fSLisandro Dalcin   Level: Advanced
86826f2ff8fSLisandro Dalcin 
86926f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
87026f2ff8fSLisandro Dalcin @*/
87126f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
87226f2ff8fSLisandro Dalcin {
87326f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
87426f2ff8fSLisandro Dalcin 
87526f2ff8fSLisandro Dalcin   PetscFunctionBegin;
87626f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
87726f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
87826f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
87926f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
88026f2ff8fSLisandro Dalcin }
88126f2ff8fSLisandro Dalcin 
88226f2ff8fSLisandro Dalcin #undef __FUNCT__
883eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
884eb284becSJed Brown /*@
885eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
886eb284becSJed Brown 
887eb284becSJed Brown   Not Collective
888eb284becSJed Brown 
889eb284becSJed Brown   Input Parameter:
890eb284becSJed Brown +  ts - timestepping context
891eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
892eb284becSJed Brown 
893eb284becSJed Brown   Options Database:
894eb284becSJed Brown .  -ts_theta_endpoint <flg>
895eb284becSJed Brown 
896eb284becSJed Brown   Level: Intermediate
897eb284becSJed Brown 
898eb284becSJed Brown .seealso: TSTHETA, TSCN
899eb284becSJed Brown @*/
900eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
901eb284becSJed Brown {
902eb284becSJed Brown   PetscErrorCode ierr;
903eb284becSJed Brown 
904eb284becSJed Brown   PetscFunctionBegin;
905eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
906eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
907eb284becSJed Brown   PetscFunctionReturn(0);
908eb284becSJed Brown }
909eb284becSJed Brown 
910f33bbcb6SJed Brown /*
911f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
912f33bbcb6SJed Brown  * The creation functions for these specializations are below.
913f33bbcb6SJed Brown  */
914f33bbcb6SJed Brown 
915f33bbcb6SJed Brown #undef __FUNCT__
916f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
917f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
918f33bbcb6SJed Brown {
919d52bd9f3SBarry Smith   PetscErrorCode ierr;
920d52bd9f3SBarry Smith 
921f33bbcb6SJed Brown   PetscFunctionBegin;
922d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
923f33bbcb6SJed Brown   PetscFunctionReturn(0);
924f33bbcb6SJed Brown }
925f33bbcb6SJed Brown 
926f33bbcb6SJed Brown /*MC
927f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
928f33bbcb6SJed Brown 
929f33bbcb6SJed Brown   Level: beginner
930f33bbcb6SJed Brown 
9314eb428fdSBarry Smith   Notes:
932c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
9334eb428fdSBarry Smith 
9344eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
9354eb428fdSBarry Smith 
936f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
937f33bbcb6SJed Brown 
938f33bbcb6SJed Brown M*/
939f33bbcb6SJed Brown #undef __FUNCT__
940f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
9418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
942f33bbcb6SJed Brown {
943f33bbcb6SJed Brown   PetscErrorCode ierr;
944f33bbcb6SJed Brown 
945f33bbcb6SJed Brown   PetscFunctionBegin;
946f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
947f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
9480c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
949f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
950f33bbcb6SJed Brown   PetscFunctionReturn(0);
951f33bbcb6SJed Brown }
952f33bbcb6SJed Brown 
953f33bbcb6SJed Brown #undef __FUNCT__
954f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
955f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
956f33bbcb6SJed Brown {
957d52bd9f3SBarry Smith   PetscErrorCode ierr;
958d52bd9f3SBarry Smith 
959f33bbcb6SJed Brown   PetscFunctionBegin;
960d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
961f33bbcb6SJed Brown   PetscFunctionReturn(0);
962f33bbcb6SJed Brown }
963f33bbcb6SJed Brown 
964f33bbcb6SJed Brown /*MC
965f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
966f33bbcb6SJed Brown 
967f33bbcb6SJed Brown   Level: beginner
968f33bbcb6SJed Brown 
969f33bbcb6SJed Brown   Notes:
9707cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
9717cf5af47SJed Brown 
9727cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
973f33bbcb6SJed Brown 
974f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
975f33bbcb6SJed Brown 
976f33bbcb6SJed Brown M*/
977f33bbcb6SJed Brown #undef __FUNCT__
978f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
9798cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
980f33bbcb6SJed Brown {
981f33bbcb6SJed Brown   PetscErrorCode ierr;
982f33bbcb6SJed Brown 
983f33bbcb6SJed Brown   PetscFunctionBegin;
984f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
985f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
986eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
9870c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
988f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
989f33bbcb6SJed Brown   PetscFunctionReturn(0);
990f33bbcb6SJed Brown }
991