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 { 10715f1b00SHong Zhang /* context for time stepping */ 111566a47fSLisandro Dalcin PetscReal stage_time; 121566a47fSLisandro Dalcin Vec X0,X,Xdot; /* Storage for stages and time derivative */ 131566a47fSLisandro Dalcin Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 141566a47fSLisandro Dalcin PetscReal Theta; 15715f1b00SHong Zhang PetscReal ptime; 16715f1b00SHong Zhang PetscReal time_step; 171566a47fSLisandro Dalcin PetscInt order; 181566a47fSLisandro Dalcin PetscBool endpoint; 191566a47fSLisandro Dalcin PetscBool extrapolate; 20715f1b00SHong Zhang TSStepStatus status; 21715f1b00SHong Zhang Vec VecCostIntegral0; /* Backup for roll-backs due to events */ 221566a47fSLisandro Dalcin 23715f1b00SHong Zhang /* context for sensitivity analysis */ 24022c081aSHong Zhang PetscInt num_tlm; /* Total number of tangent linear equations */ 25b5e0532dSHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 26b5e0532dSHong Zhang Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 27b5e0532dSHong Zhang Vec *VecsSensiTemp; /* Vector to be timed with Jacobian transpose */ 28715f1b00SHong Zhang Vec *VecsDeltaFwdSensi; /* Increment of the forward sensitivity at stage */ 29715f1b00SHong Zhang Vec VecIntegralSensiTemp; /* Working vector for forward integral sensitivity */ 30715f1b00SHong Zhang Vec VecIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 316f92202bSHong Zhang Vec *VecsFwdSensi0; /* backup for roll-backs due to events */ 326f92202bSHong Zhang Vec *VecsIntegralSensi0; /* backup for roll-backs due to events */ 336f92202bSHong Zhang Vec *VecsIntegralSensip0; /* backup for roll-backs due to events */ 341566a47fSLisandro Dalcin 35715f1b00SHong Zhang /* context for error estimation */ 361566a47fSLisandro Dalcin Vec vec_sol_prev; 371566a47fSLisandro Dalcin Vec vec_lte_work; 38316643e7SJed Brown } TS_Theta; 39316643e7SJed Brown 407445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 417445fe48SJed Brown { 427445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 437445fe48SJed Brown PetscErrorCode ierr; 447445fe48SJed Brown 457445fe48SJed Brown PetscFunctionBegin; 467445fe48SJed Brown if (X0) { 477445fe48SJed Brown if (dm && dm != ts->dm) { 480d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 497445fe48SJed Brown } else *X0 = ts->vec_sol; 507445fe48SJed Brown } 517445fe48SJed Brown if (Xdot) { 527445fe48SJed Brown if (dm && dm != ts->dm) { 530d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 547445fe48SJed Brown } else *Xdot = th->Xdot; 557445fe48SJed Brown } 567445fe48SJed Brown PetscFunctionReturn(0); 577445fe48SJed Brown } 587445fe48SJed Brown 590d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 600d0b770aSPeter Brune { 610d0b770aSPeter Brune PetscErrorCode ierr; 620d0b770aSPeter Brune 630d0b770aSPeter Brune PetscFunctionBegin; 640d0b770aSPeter Brune if (X0) { 650d0b770aSPeter Brune if (dm && dm != ts->dm) { 660d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 670d0b770aSPeter Brune } 680d0b770aSPeter Brune } 690d0b770aSPeter Brune if (Xdot) { 700d0b770aSPeter Brune if (dm && dm != ts->dm) { 710d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 720d0b770aSPeter Brune } 730d0b770aSPeter Brune } 740d0b770aSPeter Brune PetscFunctionReturn(0); 750d0b770aSPeter Brune } 760d0b770aSPeter Brune 777445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 787445fe48SJed Brown { 797445fe48SJed Brown PetscFunctionBegin; 807445fe48SJed Brown PetscFunctionReturn(0); 817445fe48SJed Brown } 827445fe48SJed Brown 837445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 847445fe48SJed Brown { 857445fe48SJed Brown TS ts = (TS)ctx; 867445fe48SJed Brown PetscErrorCode ierr; 877445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 887445fe48SJed Brown 897445fe48SJed Brown PetscFunctionBegin; 907445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 917445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 927445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 937445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 947445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 957445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 960d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 970d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 987445fe48SJed Brown PetscFunctionReturn(0); 997445fe48SJed Brown } 1007445fe48SJed Brown 101258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 102258e1594SPeter Brune { 103258e1594SPeter Brune PetscFunctionBegin; 104258e1594SPeter Brune PetscFunctionReturn(0); 105258e1594SPeter Brune } 106258e1594SPeter Brune 107258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 108258e1594SPeter Brune { 109258e1594SPeter Brune TS ts = (TS)ctx; 110258e1594SPeter Brune PetscErrorCode ierr; 111258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 112258e1594SPeter Brune 113258e1594SPeter Brune PetscFunctionBegin; 114258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 115258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 116258e1594SPeter Brune 117258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 118258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 119258e1594SPeter Brune 120258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 122258e1594SPeter Brune 123258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 124258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 125258e1594SPeter Brune PetscFunctionReturn(0); 126258e1594SPeter Brune } 127258e1594SPeter Brune 128022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 129022c081aSHong Zhang { 130022c081aSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 131022c081aSHong Zhang PetscErrorCode ierr; 132022c081aSHong Zhang 133022c081aSHong Zhang PetscFunctionBegin; 134022c081aSHong Zhang if (th->endpoint) { 135022c081aSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 136022c081aSHong Zhang if (th->Theta!=1.0) { 137022c081aSHong Zhang ierr = TSComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 138022c081aSHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 139022c081aSHong Zhang } 140022c081aSHong Zhang ierr = TSComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 141022c081aSHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 142022c081aSHong Zhang } else { 143022c081aSHong Zhang ierr = TSComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 144022c081aSHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 145022c081aSHong Zhang } 146022c081aSHong Zhang PetscFunctionReturn(0); 147022c081aSHong Zhang } 148022c081aSHong Zhang 149b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 150b1cb36f3SHong Zhang { 151b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 152b1cb36f3SHong Zhang PetscErrorCode ierr; 153b1cb36f3SHong Zhang 154b1cb36f3SHong Zhang PetscFunctionBegin; 155b1cb36f3SHong Zhang /* backup cost integral */ 156b1cb36f3SHong Zhang ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr); 157022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 158b1cb36f3SHong Zhang PetscFunctionReturn(0); 159b1cb36f3SHong Zhang } 160b1cb36f3SHong Zhang 161b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 162b1cb36f3SHong Zhang { 163b1cb36f3SHong Zhang PetscErrorCode ierr; 164b1cb36f3SHong Zhang 165b1cb36f3SHong Zhang PetscFunctionBegin; 166022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 16724655328SShri PetscFunctionReturn(0); 16824655328SShri } 16924655328SShri 1701566a47fSLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x) 1711566a47fSLisandro Dalcin { 1721566a47fSLisandro Dalcin PetscInt nits,lits; 1731566a47fSLisandro Dalcin PetscErrorCode ierr; 1741566a47fSLisandro Dalcin 1751566a47fSLisandro Dalcin PetscFunctionBegin; 1761566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1771566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1781566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1791566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1801566a47fSLisandro Dalcin PetscFunctionReturn(0); 1811566a47fSLisandro Dalcin } 1821566a47fSLisandro Dalcin 183193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 184316643e7SJed Brown { 185316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1861566a47fSLisandro Dalcin PetscInt rejections = 0; 1874957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 1881566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 189051f2191SLisandro Dalcin PetscErrorCode ierr; 190316643e7SJed Brown 191316643e7SJed Brown PetscFunctionBegin; 1921566a47fSLisandro Dalcin if (!ts->steprollback) { 1932ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 1943b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 1951566a47fSLisandro Dalcin } 1961566a47fSLisandro Dalcin 1971566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 1981566a47fSLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 199715f1b00SHong Zhang 2001566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 2011566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 202316643e7SJed Brown 203be5899b3SLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 204fecfb714SLisandro Dalcin if (th->extrapolate && !ts->steprestart) { 205be5899b3SLisandro Dalcin ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 206be5899b3SLisandro Dalcin } 207eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 208eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2091566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2101566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2111566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2121566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2131566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 214eb284becSJed Brown } 215be5899b3SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 2161566a47fSLisandro Dalcin ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 217be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 2181566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 219fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 220051f2191SLisandro Dalcin 221051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2221566a47fSLisandro Dalcin if (th->endpoint) { 2231566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2241566a47fSLisandro Dalcin } else { 225cb3a5882SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 2261566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2271566a47fSLisandro Dalcin } 2281566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2291566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2301566a47fSLisandro Dalcin if (!accept) { 2311566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 232be5899b3SLisandro Dalcin ts->time_step = next_time_step; 233be5899b3SLisandro Dalcin goto reject_step; 234051f2191SLisandro Dalcin } 2353b1890cdSShri Abhyankar 236715f1b00SHong Zhang if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 237b1cb36f3SHong Zhang th->ptime = ts->ptime; 238b1cb36f3SHong Zhang th->time_step = ts->time_step; 23917145e75SHong Zhang } 2401566a47fSLisandro Dalcin 2412b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 242cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 243051f2191SLisandro Dalcin break; 244051f2191SLisandro Dalcin 245051f2191SLisandro Dalcin reject_step: 246fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 2471566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 248051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2491566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2503b1890cdSShri Abhyankar } 2513b1890cdSShri Abhyankar } 252316643e7SJed Brown PetscFunctionReturn(0); 253316643e7SJed Brown } 254316643e7SJed Brown 25542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 2562ca6e920SHong Zhang { 2572ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 258b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 2592ca6e920SHong Zhang PetscInt nadj; 2602ca6e920SHong Zhang PetscErrorCode ierr; 2612ca6e920SHong Zhang Mat J,Jp; 2622ca6e920SHong Zhang KSP ksp; 2632ca6e920SHong Zhang PetscReal shift; 2642ca6e920SHong Zhang 2652ca6e920SHong Zhang PetscFunctionBegin; 2662ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 267302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 2682ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 269b5e0532dSHong Zhang 270b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 271022c081aSHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 272b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 273022c081aSHong Zhang th->time_step = -ts->time_step; 2742ca6e920SHong Zhang 275a4cab896SHong Zhang /* Build RHS */ 2762c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 27705755b9cSHong Zhang if (th->endpoint) { 278022c081aSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,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); 285022c081aSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->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 */ 292022c081aSHong Zhang shift = 1./(th->Theta*th->time_step); 2933c54f92cSHong Zhang if (th->endpoint) { 294022c081aSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,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.) { 308022c081aSHong Zhang shift = 1./((th->Theta-1.)*th->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); 3128263dbbfSHong Zhang } 313abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 314b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 3152c39e106SBarry Smith if (ts->vec_costintegral) { 31636eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 31736eaed60SHong Zhang } 3183c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 3192ca6e920SHong Zhang } 320b5e0532dSHong Zhang } else { /* backward Euler */ 321b5e0532dSHong Zhang shift = 0.0; 322022c081aSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 323b5e0532dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 324b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 325022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 326b5e0532dSHong Zhang if (ts->vec_costintegral) { 327022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 328b5e0532dSHong Zhang } 329b5e0532dSHong Zhang } 330b5e0532dSHong Zhang } 3313fd52205SHong Zhang 3323fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 333022c081aSHong Zhang ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 334abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 335b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 336022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3373fd52205SHong Zhang } 338b5e0532dSHong Zhang if (th->Theta!=1.) { 339b5e0532dSHong Zhang ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 340abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 341b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 342022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 343b5e0532dSHong Zhang } 3443fd52205SHong Zhang } 3452c39e106SBarry Smith if (ts->vec_costintegral) { 346022c081aSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 347abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 348022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 34936eaed60SHong Zhang } 350b5e0532dSHong Zhang if (th->Theta!=1.) { 351b5e0532dSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 352abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 353022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 35436eaed60SHong Zhang } 35536eaed60SHong Zhang } 3563fd52205SHong Zhang } 357b5e0532dSHong Zhang } 3583fd52205SHong Zhang } else { /* one-stage case */ 3593c54f92cSHong Zhang shift = 0.0; 360a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 3612c39e106SBarry Smith if (ts->vec_costintegral) { 362d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 3633c54f92cSHong Zhang } 364abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 365b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 366022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 3672c39e106SBarry Smith if (ts->vec_costintegral) { 368022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 36936eaed60SHong Zhang } 3702ca6e920SHong Zhang } 3713fd52205SHong Zhang if (ts->vecs_sensip) { 3725bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 373abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 374b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 375022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3763fd52205SHong Zhang } 3772c39e106SBarry Smith if (ts->vec_costintegral) { 378d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 379abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 380022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 38136eaed60SHong Zhang } 38236eaed60SHong Zhang } 3833fd52205SHong Zhang } 3842ca6e920SHong Zhang } 3852ca6e920SHong Zhang 3862ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 3872ca6e920SHong Zhang PetscFunctionReturn(0); 3882ca6e920SHong Zhang } 3892ca6e920SHong Zhang 390cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 391cd652676SJed Brown { 392cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 393be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 394cd652676SJed Brown PetscErrorCode ierr; 395cd652676SJed Brown 396cd652676SJed Brown PetscFunctionBegin; 397a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 398be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 399be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 400cd652676SJed Brown PetscFunctionReturn(0); 401cd652676SJed Brown } 402cd652676SJed Brown 4031566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 4041566a47fSLisandro Dalcin { 4051566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 4061566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 4071566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 4087453f775SEmil Constantinescu PetscReal wltea,wlter; 4091566a47fSLisandro Dalcin PetscErrorCode ierr; 4101566a47fSLisandro Dalcin 4111566a47fSLisandro Dalcin PetscFunctionBegin; 4122ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 4131566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 414fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 4151566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 416fecfb714SLisandro Dalcin { 417be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 418be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 4191566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 4201566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 4211566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 4221566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 4231566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 4247453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 4251566a47fSLisandro Dalcin } 4261566a47fSLisandro Dalcin if (order) *order = 2; 4271566a47fSLisandro Dalcin PetscFunctionReturn(0); 4281566a47fSLisandro Dalcin } 4291566a47fSLisandro Dalcin 4301566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 4311566a47fSLisandro Dalcin { 4321566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 433022c081aSHong Zhang PetscInt ntlm,ncost; 4341566a47fSLisandro Dalcin PetscErrorCode ierr; 4351566a47fSLisandro Dalcin 4361566a47fSLisandro Dalcin PetscFunctionBegin; 4371566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 4381566a47fSLisandro Dalcin if (ts->vec_costintegral && ts->costintegralfwd) { 4391566a47fSLisandro Dalcin ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 4401566a47fSLisandro Dalcin } 441715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 4426f92202bSHong Zhang 443022c081aSHong Zhang for (ntlm=0;ntlm<th->num_tlm;ntlm++) { 444022c081aSHong Zhang ierr = VecCopy(th->VecsFwdSensi0[ntlm],ts->vecs_fwdsensipacked[ntlm]);CHKERRQ(ierr); 4456f92202bSHong Zhang } 4466f92202bSHong Zhang if (ts->vecs_integral_sensi) { 4476f92202bSHong Zhang for (ncost=0;ncost<ts->numcost;ncost++) { 4486f92202bSHong Zhang ierr = VecCopy(th->VecsIntegralSensi0[ncost],ts->vecs_integral_sensi[ncost]);CHKERRQ(ierr); 4496f92202bSHong Zhang } 4506f92202bSHong Zhang } 4516f92202bSHong Zhang if (ts->vecs_integral_sensip) { 4526f92202bSHong Zhang for (ncost=0;ncost<ts->numcost;ncost++) { 4536f92202bSHong Zhang ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 4546f92202bSHong Zhang } 4556f92202bSHong Zhang } 456715f1b00SHong Zhang PetscFunctionReturn(0); 457715f1b00SHong Zhang } 458715f1b00SHong Zhang 459022c081aSHong Zhang static PetscErrorCode TSThetaIntegrandDerivative(TS ts,PetscInt numtlm,Vec DRncostDY,Vec* s,Vec DRncostDP,Vec VecIntegrandDerivative) 460715f1b00SHong Zhang { 461022c081aSHong Zhang PetscInt ntlm,low,high; 462715f1b00SHong Zhang PetscScalar *v; 463715f1b00SHong Zhang PetscErrorCode ierr; 464715f1b00SHong Zhang 465715f1b00SHong Zhang PetscFunctionBegin; 466715f1b00SHong Zhang 467022c081aSHong Zhang ierr = VecGetOwnershipRange(VecIntegrandDerivative,&low,&high);CHKERRQ(ierr); 468715f1b00SHong Zhang ierr = VecGetArray(VecIntegrandDerivative,&v);CHKERRQ(ierr); 469022c081aSHong Zhang for (ntlm=low; ntlm<high; ntlm++) { 470022c081aSHong Zhang ierr = VecDot(DRncostDY,s[ntlm],&v[ntlm-low]);CHKERRQ(ierr); 471715f1b00SHong Zhang } 472715f1b00SHong Zhang ierr = VecRestoreArray(VecIntegrandDerivative,&v);CHKERRQ(ierr); 473715f1b00SHong Zhang if (DRncostDP) { 474715f1b00SHong Zhang ierr = VecAXPY(VecIntegrandDerivative,1.,DRncostDP);CHKERRQ(ierr); 475715f1b00SHong Zhang } 476715f1b00SHong Zhang PetscFunctionReturn(0); 477715f1b00SHong Zhang } 478715f1b00SHong Zhang 479715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 480715f1b00SHong Zhang { 481715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 482715f1b00SHong Zhang Vec *VecsDeltaFwdSensi = th->VecsDeltaFwdSensi; 483022c081aSHong Zhang PetscInt ntlm,ncost,np; 484715f1b00SHong Zhang KSP ksp; 485715f1b00SHong Zhang Mat J,Jp; 486715f1b00SHong Zhang PetscReal shift; 487715f1b00SHong Zhang PetscErrorCode ierr; 488715f1b00SHong Zhang 489715f1b00SHong Zhang PetscFunctionBegin; 490022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 491022c081aSHong Zhang ierr = VecCopy(ts->vecs_fwdsensipacked[ntlm],th->VecsFwdSensi0[ntlm]);CHKERRQ(ierr); 4926f92202bSHong Zhang } 4936f92202bSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 494022c081aSHong Zhang if (ts->vecs_integral_sensi) { 4956f92202bSHong Zhang ierr = VecCopy(ts->vecs_integral_sensi[ncost],th->VecsIntegralSensi0[ncost]);CHKERRQ(ierr); 4966f92202bSHong Zhang } 4976f92202bSHong Zhang if (ts->vecs_integral_sensip) { 4986f92202bSHong Zhang ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 4996f92202bSHong Zhang } 5006f92202bSHong Zhang } 5016f92202bSHong Zhang 502715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 503715f1b00SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 504715f1b00SHong Zhang 505715f1b00SHong Zhang /* Build RHS */ 506715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 507715f1b00SHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 508715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 509715f1b00SHong Zhang 510022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 5117f79407eSBarry Smith ierr = MatMult(J,ts->vecs_fwdsensipacked[ntlm],VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr); 5127f79407eSBarry Smith ierr = VecScale(VecsDeltaFwdSensi[ntlm],(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 513022c081aSHong Zhang } 514022c081aSHong Zhang /* Add the f_p forcing terms */ 515715f1b00SHong Zhang ierr = TSForwardComputeRHSJacobianP(ts,th->ptime,th->X0,ts->vecs_jacp);CHKERRQ(ierr); 516022c081aSHong Zhang for (np=0; np<ts->num_parameters; np++) { 517022c081aSHong Zhang ierr = VecAXPY(VecsDeltaFwdSensi[np+ts->num_initialvalues],(1.-th->Theta)/th->Theta,ts->vecs_jacp[np]);CHKERRQ(ierr); 518715f1b00SHong Zhang } 519715f1b00SHong Zhang ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->vecs_jacp);CHKERRQ(ierr); 520022c081aSHong Zhang for (np=0; np<ts->num_parameters; np++) { 521022c081aSHong Zhang ierr = VecAXPY(VecsDeltaFwdSensi[np+ts->num_initialvalues],1,ts->vecs_jacp[np]);CHKERRQ(ierr); 522715f1b00SHong Zhang } 523715f1b00SHong Zhang } else { /* 1-stage method */ 524715f1b00SHong Zhang shift = 0.0; 525715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 526022c081aSHong Zhang for (ntlm=0; ntlm<ts->num_parameters; ntlm++) { 5277f79407eSBarry Smith ierr = MatMult(J,ts->vecs_fwdsensipacked[ntlm],VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr); 5287f79407eSBarry Smith ierr = VecScale(VecsDeltaFwdSensi[ntlm],-1);CHKERRQ(ierr); 529715f1b00SHong Zhang } 530022c081aSHong Zhang /* Add the f_p forcing terms */ 531715f1b00SHong Zhang ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,th->X,ts->vecs_jacp);CHKERRQ(ierr); 532022c081aSHong Zhang for (np=0; np<ts->num_parameters; np++) { 533022c081aSHong Zhang ierr = VecAXPY(VecsDeltaFwdSensi[np+ts->num_initialvalues],1,ts->vecs_jacp[np]);CHKERRQ(ierr); 534715f1b00SHong Zhang } 535715f1b00SHong Zhang } 536715f1b00SHong Zhang 537715f1b00SHong Zhang /* Build LHS */ 538715f1b00SHong Zhang shift = 1/(th->Theta*th->time_step); 539715f1b00SHong Zhang if (th->endpoint) { 540715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 541715f1b00SHong Zhang } else { 542715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 543715f1b00SHong Zhang } 544715f1b00SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 545715f1b00SHong Zhang 546715f1b00SHong Zhang /* 547715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 548715f1b00SHong Zhang drdy|t_n*S(t_n) + drdp|t_n 549715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 550715f1b00SHong Zhang It is important to note that sensitivitesi to parameters (sensip) and sensitivities to initial values (sensi) are independent of each other, but integral sensip relies on sensip and integral sensi requires sensi 551715f1b00SHong Zhang */ 552715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 553715f1b00SHong Zhang if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) { 554715f1b00SHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr); 555715f1b00SHong Zhang } 556715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 557715f1b00SHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 558715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 5597f79407eSBarry Smith ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],&ts->vecs_fwdsensipacked[ts->num_initialvalues],ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 560715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 561715f1b00SHong Zhang } 562715f1b00SHong Zhang } 563715f1b00SHong Zhang if (ts->vecs_integral_sensi) { 564715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 5657f79407eSBarry Smith ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensipacked,NULL,th->VecIntegralSensiTemp);CHKERRQ(ierr); 566715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensiTemp);CHKERRQ(ierr); 567715f1b00SHong Zhang } 568715f1b00SHong Zhang } 569715f1b00SHong Zhang } 570715f1b00SHong Zhang 571715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 572022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 573715f1b00SHong Zhang if (th->endpoint) { 5747f79407eSBarry Smith ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ntlm],ts->vecs_fwdsensipacked[ntlm]);CHKERRQ(ierr); 575715f1b00SHong Zhang } else { 5767f79407eSBarry Smith ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ntlm],VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr); 5777f79407eSBarry Smith ierr = VecAXPY(ts->vecs_fwdsensipacked[ntlm],1.,VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr); 578715f1b00SHong Zhang } 579715f1b00SHong Zhang } 580715f1b00SHong Zhang /* Evaluate the second stage of integral gradients with the 2-stage method: 581715f1b00SHong Zhang drdy|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 582715f1b00SHong Zhang */ 583715f1b00SHong Zhang if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) { 584715f1b00SHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 585715f1b00SHong Zhang } 586715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 587715f1b00SHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 588715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 5897f79407eSBarry Smith ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],&ts->vecs_fwdsensipacked[ts->num_initialvalues],ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 590715f1b00SHong Zhang if (th->endpoint) { 591715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 592715f1b00SHong Zhang } else { 593715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 594715f1b00SHong Zhang } 595715f1b00SHong Zhang } 596715f1b00SHong Zhang } 597715f1b00SHong Zhang if (ts->vecs_integral_sensi) { 598715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 5997f79407eSBarry Smith ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensipacked,NULL,th->VecIntegralSensiTemp);CHKERRQ(ierr); 600715f1b00SHong Zhang if (th->endpoint) { 601715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*th->Theta,th->VecIntegralSensiTemp);CHKERRQ(ierr); 602715f1b00SHong Zhang } else { 603715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step,th->VecIntegralSensiTemp);CHKERRQ(ierr); 604715f1b00SHong Zhang } 605715f1b00SHong Zhang } 606715f1b00SHong Zhang } 607715f1b00SHong Zhang 608715f1b00SHong Zhang if (!th->endpoint) { /* VecsDeltaFwdSensip are the increment for the stage values for the 1-stage method */ 609022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 6107f79407eSBarry Smith ierr = VecAXPY(ts->vecs_fwdsensipacked[ntlm],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr); 611715f1b00SHong Zhang } 612715f1b00SHong Zhang } 6131566a47fSLisandro Dalcin PetscFunctionReturn(0); 6141566a47fSLisandro Dalcin } 6151566a47fSLisandro Dalcin 616316643e7SJed Brown /*------------------------------------------------------------*/ 617277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 618316643e7SJed Brown { 619316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 620316643e7SJed Brown PetscErrorCode ierr; 621316643e7SJed Brown 622316643e7SJed Brown PetscFunctionBegin; 6236bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 6246bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 6253b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 626eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 6271566a47fSLisandro Dalcin 6281566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 6291566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 6301566a47fSLisandro Dalcin 631b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 632715f1b00SHong Zhang if (ts->forward_solve) { 633022c081aSHong Zhang ierr = VecDestroyVecs(th->num_tlm,&th->VecsDeltaFwdSensi);CHKERRQ(ierr); 634022c081aSHong Zhang ierr = VecDestroyVecs(th->num_tlm,&th->VecsFwdSensi0);CHKERRQ(ierr); 635715f1b00SHong Zhang if (ts->vecs_integral_sensi) { 636715f1b00SHong Zhang ierr = VecDestroy(&th->VecIntegralSensiTemp);CHKERRQ(ierr); 6376f92202bSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensi0);CHKERRQ(ierr); 638715f1b00SHong Zhang } 639715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 640715f1b00SHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 6416f92202bSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 642715f1b00SHong Zhang } 643715f1b00SHong Zhang } 644b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 645b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 646b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 647715f1b00SHong Zhang 648277b19d0SLisandro Dalcin PetscFunctionReturn(0); 649277b19d0SLisandro Dalcin } 650277b19d0SLisandro Dalcin 651277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 652277b19d0SLisandro Dalcin { 653277b19d0SLisandro Dalcin PetscErrorCode ierr; 654277b19d0SLisandro Dalcin 655277b19d0SLisandro Dalcin PetscFunctionBegin; 656277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 657*b3a6b972SJed Brown if (ts->dm) { 658*b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 659*b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 660*b3a6b972SJed Brown } 661277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 662bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 663bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 664bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 665bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 666316643e7SJed Brown PetscFunctionReturn(0); 667316643e7SJed Brown } 668316643e7SJed Brown 669316643e7SJed Brown /* 670316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 6712b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 672316643e7SJed Brown */ 6730f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 674316643e7SJed Brown { 675316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 676316643e7SJed Brown PetscErrorCode ierr; 6777445fe48SJed Brown Vec X0,Xdot; 6787445fe48SJed Brown DM dm,dmsave; 6791566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 680316643e7SJed Brown 681316643e7SJed Brown PetscFunctionBegin; 6827445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 6835a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 6847445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 685b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 6867445fe48SJed Brown 6877445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 6887445fe48SJed Brown dmsave = ts->dm; 6897445fe48SJed Brown ts->dm = dm; 6907445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 6917445fe48SJed Brown ts->dm = dmsave; 6920d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 693316643e7SJed Brown PetscFunctionReturn(0); 694316643e7SJed Brown } 695316643e7SJed Brown 696d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 697316643e7SJed Brown { 698316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 699316643e7SJed Brown PetscErrorCode ierr; 7007445fe48SJed Brown Vec Xdot; 7017445fe48SJed Brown DM dm,dmsave; 7021566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 703316643e7SJed Brown 704316643e7SJed Brown PetscFunctionBegin; 7057445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 706be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 7070298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 7087445fe48SJed Brown 7097445fe48SJed Brown dmsave = ts->dm; 7107445fe48SJed Brown ts->dm = dm; 711d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 7127445fe48SJed Brown ts->dm = dmsave; 7130298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 714316643e7SJed Brown PetscFunctionReturn(0); 715316643e7SJed Brown } 716316643e7SJed Brown 717715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 718715f1b00SHong Zhang { 719715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 720715f1b00SHong Zhang PetscErrorCode ierr; 721715f1b00SHong Zhang 722715f1b00SHong Zhang PetscFunctionBegin; 723022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 724022c081aSHong Zhang th->num_tlm = ts->num_parameters + ts->num_initialvalues; 725022c081aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_fwdsensipacked[0],th->num_tlm,&th->VecsDeltaFwdSensi);CHKERRQ(ierr); 726715f1b00SHong Zhang if (ts->vecs_integral_sensi) { 727715f1b00SHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensi[0],&th->VecIntegralSensiTemp);CHKERRQ(ierr); 728715f1b00SHong Zhang } 729715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 730715f1b00SHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 731715f1b00SHong Zhang } 7326f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 733022c081aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_fwdsensipacked[0],th->num_tlm,&th->VecsFwdSensi0);CHKERRQ(ierr); 7346f92202bSHong Zhang if (ts->vecs_integral_sensi) { 7356f92202bSHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensi[0],ts->numcost,&th->VecsIntegralSensi0);CHKERRQ(ierr); 7366f92202bSHong Zhang } 7376f92202bSHong Zhang if (ts->vecs_integral_sensip) { 7386f92202bSHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 7396f92202bSHong Zhang } 740715f1b00SHong Zhang PetscFunctionReturn(0); 741715f1b00SHong Zhang } 742715f1b00SHong Zhang 743316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 744316643e7SJed Brown { 745316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 7462ffb9264SLisandro Dalcin PetscBool match; 747316643e7SJed Brown PetscErrorCode ierr; 748316643e7SJed Brown 749316643e7SJed Brown PetscFunctionBegin; 750b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 751b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 752b1cb36f3SHong Zhang } 75339d13666SHong Zhang if (!th->X) { 754316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 75539d13666SHong Zhang } 75639d13666SHong Zhang if (!th->Xdot) { 757316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 75839d13666SHong Zhang } 75939d13666SHong Zhang if (!th->X0) { 7603b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 76139d13666SHong Zhang } 7621566a47fSLisandro Dalcin if (th->endpoint) { 7631566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 7647445fe48SJed Brown } 7653b1890cdSShri Abhyankar 7661566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 7671566a47fSLisandro Dalcin 7681566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 7691566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 7701566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 7711566a47fSLisandro Dalcin 7721566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 7731566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 7742ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 7752ffb9264SLisandro Dalcin if (!match) { 7761566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 7771566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 7783b1890cdSShri Abhyankar } 7791566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 780b4dd3bf9SBarry Smith PetscFunctionReturn(0); 781b4dd3bf9SBarry Smith } 7820c7376e5SHong Zhang 783b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 784b4dd3bf9SBarry Smith 785b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 786b4dd3bf9SBarry Smith { 787b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 788b4dd3bf9SBarry Smith PetscErrorCode ierr; 789b4dd3bf9SBarry Smith 790b4dd3bf9SBarry Smith PetscFunctionBegin; 791b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 7922ca6e920SHong Zhang if(ts->vecs_sensip) { 793b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 7942ca6e920SHong Zhang } 795b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 796316643e7SJed Brown PetscFunctionReturn(0); 797316643e7SJed Brown } 798316643e7SJed Brown 7994416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 800316643e7SJed Brown { 801316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 802316643e7SJed Brown PetscErrorCode ierr; 803316643e7SJed Brown 804316643e7SJed Brown PetscFunctionBegin; 805e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 806316643e7SJed Brown { 8070298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 8080298fd71SBarry 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); 80903842d09SLisandro 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); 810316643e7SJed Brown } 811316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 812316643e7SJed Brown PetscFunctionReturn(0); 813316643e7SJed Brown } 814316643e7SJed Brown 815316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 816316643e7SJed Brown { 817316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 818ace3abfcSBarry Smith PetscBool iascii; 819316643e7SJed Brown PetscErrorCode ierr; 820316643e7SJed Brown 821316643e7SJed Brown PetscFunctionBegin; 822251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 823316643e7SJed Brown if (iascii) { 8247c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 825ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 826316643e7SJed Brown } 827316643e7SJed Brown PetscFunctionReturn(0); 828316643e7SJed Brown } 829316643e7SJed Brown 830560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 8310de4c49aSJed Brown { 8320de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 8330de4c49aSJed Brown 8340de4c49aSJed Brown PetscFunctionBegin; 8350de4c49aSJed Brown *theta = th->Theta; 8360de4c49aSJed Brown PetscFunctionReturn(0); 8370de4c49aSJed Brown } 8380de4c49aSJed Brown 839560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 8400de4c49aSJed Brown { 8410de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 8420de4c49aSJed Brown 8430de4c49aSJed Brown PetscFunctionBegin; 8447c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 8450de4c49aSJed Brown th->Theta = theta; 8461566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 8470de4c49aSJed Brown PetscFunctionReturn(0); 8480de4c49aSJed Brown } 849eb284becSJed Brown 850560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 85126f2ff8fSLisandro Dalcin { 85226f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 85326f2ff8fSLisandro Dalcin 85426f2ff8fSLisandro Dalcin PetscFunctionBegin; 85526f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 85626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 85726f2ff8fSLisandro Dalcin } 85826f2ff8fSLisandro Dalcin 859560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 860eb284becSJed Brown { 861eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 862eb284becSJed Brown 863eb284becSJed Brown PetscFunctionBegin; 864eb284becSJed Brown th->endpoint = flg; 865eb284becSJed Brown PetscFunctionReturn(0); 866eb284becSJed Brown } 8670de4c49aSJed Brown 868f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 869f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 870f9c1d6abSBarry Smith { 871f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 872f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 8733fd8ae06SJed Brown const PetscReal one = 1.0; 874f9c1d6abSBarry Smith 875f9c1d6abSBarry Smith PetscFunctionBegin; 8763fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 877f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 878f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 879f9c1d6abSBarry Smith PetscFunctionReturn(0); 880f9c1d6abSBarry Smith } 881f9c1d6abSBarry Smith #endif 882f9c1d6abSBarry Smith 88342682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 88442682096SHong Zhang { 88542682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 88642682096SHong Zhang 88742682096SHong Zhang PetscFunctionBegin; 8881566a47fSLisandro Dalcin if (ns) *ns = 1; 8891566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 89042682096SHong Zhang PetscFunctionReturn(0); 89142682096SHong Zhang } 892f9c1d6abSBarry Smith 893316643e7SJed Brown /* ------------------------------------------------------------ */ 894316643e7SJed Brown /*MC 89596f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 896316643e7SJed Brown 897316643e7SJed Brown Level: beginner 898316643e7SJed Brown 8994eb428fdSBarry Smith Options Database: 900c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 901c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 90203842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 9034eb428fdSBarry Smith 904eb284becSJed Brown Notes: 9050c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 9060c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 9074eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 9084eb428fdSBarry Smith 909eb284becSJed Brown This method can be applied to DAE. 910eb284becSJed Brown 911eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 912eb284becSJed Brown 913eb284becSJed Brown .vb 914eb284becSJed Brown Theta | Theta 915eb284becSJed Brown ------------- 916eb284becSJed Brown | 1 917eb284becSJed Brown .ve 918eb284becSJed Brown 919eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 920eb284becSJed Brown 921eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 922eb284becSJed Brown 923eb284becSJed Brown .vb 924eb284becSJed Brown 0 | 0 0 925eb284becSJed Brown 1 | 1-Theta Theta 926eb284becSJed Brown ------------------- 927eb284becSJed Brown | 1-Theta Theta 928eb284becSJed Brown .ve 929eb284becSJed Brown 930eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 931eb284becSJed Brown 932eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 933eb284becSJed Brown 934eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 935eb284becSJed Brown 9364eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 937eb284becSJed Brown 938eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 939316643e7SJed Brown 940316643e7SJed Brown M*/ 9418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 942316643e7SJed Brown { 943316643e7SJed Brown TS_Theta *th; 944316643e7SJed Brown PetscErrorCode ierr; 945316643e7SJed Brown 946316643e7SJed Brown PetscFunctionBegin; 947277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 948316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 949316643e7SJed Brown ts->ops->view = TSView_Theta; 950316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 95142f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 952316643e7SJed Brown ts->ops->step = TSStep_Theta; 953cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 9541566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 95524655328SShri ts->ops->rollback = TSRollBack_Theta; 956316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 9570f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 9580f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 959f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 960f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 961f9c1d6abSBarry Smith #endif 96242682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 96342f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 964b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 965b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 9662ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 967715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 968715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 969316643e7SJed Brown 970efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 971efd4aadfSBarry Smith 972b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 973316643e7SJed Brown ts->data = (void*)th; 974316643e7SJed Brown 975715f1b00SHong Zhang th->VecsDeltaLam = NULL; 976715f1b00SHong Zhang th->VecsDeltaMu = NULL; 977715f1b00SHong Zhang th->VecsSensiTemp = NULL; 978715f1b00SHong Zhang 9796f700aefSJed Brown th->extrapolate = PETSC_FALSE; 980316643e7SJed Brown th->Theta = 0.5; 9811566a47fSLisandro Dalcin th->order = 2; 982bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 983bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 984bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 985bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 986316643e7SJed Brown PetscFunctionReturn(0); 987316643e7SJed Brown } 9880de4c49aSJed Brown 9890de4c49aSJed Brown /*@ 9900de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 9910de4c49aSJed Brown 9920de4c49aSJed Brown Not Collective 9930de4c49aSJed Brown 9940de4c49aSJed Brown Input Parameter: 9950de4c49aSJed Brown . ts - timestepping context 9960de4c49aSJed Brown 9970de4c49aSJed Brown Output Parameter: 9980de4c49aSJed Brown . theta - stage abscissa 9990de4c49aSJed Brown 10000de4c49aSJed Brown Note: 10010de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 10020de4c49aSJed Brown 10030de4c49aSJed Brown Level: Advanced 10040de4c49aSJed Brown 10050de4c49aSJed Brown .seealso: TSThetaSetTheta() 10060de4c49aSJed Brown @*/ 10077087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 10080de4c49aSJed Brown { 10094ac538c5SBarry Smith PetscErrorCode ierr; 10100de4c49aSJed Brown 10110de4c49aSJed Brown PetscFunctionBegin; 1012afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10130de4c49aSJed Brown PetscValidPointer(theta,2); 10144ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 10150de4c49aSJed Brown PetscFunctionReturn(0); 10160de4c49aSJed Brown } 10170de4c49aSJed Brown 10180de4c49aSJed Brown /*@ 10190de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 10200de4c49aSJed Brown 10210de4c49aSJed Brown Not Collective 10220de4c49aSJed Brown 10230de4c49aSJed Brown Input Parameter: 10240de4c49aSJed Brown + ts - timestepping context 10250de4c49aSJed Brown - theta - stage abscissa 10260de4c49aSJed Brown 10270de4c49aSJed Brown Options Database: 10280de4c49aSJed Brown . -ts_theta_theta <theta> 10290de4c49aSJed Brown 10300de4c49aSJed Brown Level: Intermediate 10310de4c49aSJed Brown 10320de4c49aSJed Brown .seealso: TSThetaGetTheta() 10330de4c49aSJed Brown @*/ 10347087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 10350de4c49aSJed Brown { 10364ac538c5SBarry Smith PetscErrorCode ierr; 10370de4c49aSJed Brown 10380de4c49aSJed Brown PetscFunctionBegin; 1039afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10404ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 10410de4c49aSJed Brown PetscFunctionReturn(0); 10420de4c49aSJed Brown } 1043f33bbcb6SJed Brown 104426f2ff8fSLisandro Dalcin /*@ 104526f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 104626f2ff8fSLisandro Dalcin 104726f2ff8fSLisandro Dalcin Not Collective 104826f2ff8fSLisandro Dalcin 104926f2ff8fSLisandro Dalcin Input Parameter: 105026f2ff8fSLisandro Dalcin . ts - timestepping context 105126f2ff8fSLisandro Dalcin 105226f2ff8fSLisandro Dalcin Output Parameter: 105326f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 105426f2ff8fSLisandro Dalcin 105526f2ff8fSLisandro Dalcin Level: Advanced 105626f2ff8fSLisandro Dalcin 105726f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 105826f2ff8fSLisandro Dalcin @*/ 105926f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 106026f2ff8fSLisandro Dalcin { 106126f2ff8fSLisandro Dalcin PetscErrorCode ierr; 106226f2ff8fSLisandro Dalcin 106326f2ff8fSLisandro Dalcin PetscFunctionBegin; 106426f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 106526f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1066163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 106726f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 106826f2ff8fSLisandro Dalcin } 106926f2ff8fSLisandro Dalcin 1070eb284becSJed Brown /*@ 1071eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1072eb284becSJed Brown 1073eb284becSJed Brown Not Collective 1074eb284becSJed Brown 1075eb284becSJed Brown Input Parameter: 1076eb284becSJed Brown + ts - timestepping context 1077eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1078eb284becSJed Brown 1079eb284becSJed Brown Options Database: 1080eb284becSJed Brown . -ts_theta_endpoint <flg> 1081eb284becSJed Brown 1082eb284becSJed Brown Level: Intermediate 1083eb284becSJed Brown 1084eb284becSJed Brown .seealso: TSTHETA, TSCN 1085eb284becSJed Brown @*/ 1086eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1087eb284becSJed Brown { 1088eb284becSJed Brown PetscErrorCode ierr; 1089eb284becSJed Brown 1090eb284becSJed Brown PetscFunctionBegin; 1091eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1092eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1093eb284becSJed Brown PetscFunctionReturn(0); 1094eb284becSJed Brown } 1095eb284becSJed Brown 1096f33bbcb6SJed Brown /* 1097f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1098f33bbcb6SJed Brown * The creation functions for these specializations are below. 1099f33bbcb6SJed Brown */ 1100f33bbcb6SJed Brown 11011566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 11021566a47fSLisandro Dalcin { 11031566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 11041566a47fSLisandro Dalcin PetscErrorCode ierr; 11051566a47fSLisandro Dalcin 11061566a47fSLisandro Dalcin PetscFunctionBegin; 11071566a47fSLisandro 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"); 11081566a47fSLisandro 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"); 11091566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 11101566a47fSLisandro Dalcin PetscFunctionReturn(0); 11111566a47fSLisandro Dalcin } 11121566a47fSLisandro Dalcin 1113f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1114f33bbcb6SJed Brown { 1115f33bbcb6SJed Brown PetscFunctionBegin; 1116f33bbcb6SJed Brown PetscFunctionReturn(0); 1117f33bbcb6SJed Brown } 1118f33bbcb6SJed Brown 1119f33bbcb6SJed Brown /*MC 1120f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1121f33bbcb6SJed Brown 1122f33bbcb6SJed Brown Level: beginner 1123f33bbcb6SJed Brown 11244eb428fdSBarry Smith Notes: 1125c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 11264eb428fdSBarry Smith 11271566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 11284eb428fdSBarry Smith 1129f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1130f33bbcb6SJed Brown 1131f33bbcb6SJed Brown M*/ 11328cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1133f33bbcb6SJed Brown { 1134f33bbcb6SJed Brown PetscErrorCode ierr; 1135f33bbcb6SJed Brown 1136f33bbcb6SJed Brown PetscFunctionBegin; 1137f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1138f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 11391566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 11400c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1141f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1142f33bbcb6SJed Brown PetscFunctionReturn(0); 1143f33bbcb6SJed Brown } 1144f33bbcb6SJed Brown 11451566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 11461566a47fSLisandro Dalcin { 11471566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 11481566a47fSLisandro Dalcin PetscErrorCode ierr; 11491566a47fSLisandro Dalcin 11501566a47fSLisandro Dalcin PetscFunctionBegin; 11511566a47fSLisandro 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"); 11521566a47fSLisandro 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"); 11531566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 11541566a47fSLisandro Dalcin PetscFunctionReturn(0); 11551566a47fSLisandro Dalcin } 11561566a47fSLisandro Dalcin 1157f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1158f33bbcb6SJed Brown { 1159f33bbcb6SJed Brown PetscFunctionBegin; 1160f33bbcb6SJed Brown PetscFunctionReturn(0); 1161f33bbcb6SJed Brown } 1162f33bbcb6SJed Brown 1163f33bbcb6SJed Brown /*MC 1164f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1165f33bbcb6SJed Brown 1166f33bbcb6SJed Brown Level: beginner 1167f33bbcb6SJed Brown 1168f33bbcb6SJed Brown Notes: 11697cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 11707cf5af47SJed Brown 11717cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1172f33bbcb6SJed Brown 1173f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1174f33bbcb6SJed Brown 1175f33bbcb6SJed Brown M*/ 11768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1177f33bbcb6SJed Brown { 1178f33bbcb6SJed Brown PetscErrorCode ierr; 1179f33bbcb6SJed Brown 1180f33bbcb6SJed Brown PetscFunctionBegin; 1181f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1182f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1183eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 11840c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1185f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1186f33bbcb6SJed Brown PetscFunctionReturn(0); 1187f33bbcb6SJed Brown } 1188