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 */ 2713b1b0ffSHong Zhang Vec *VecsSensiTemp; /* Vector to be multiplied with Jacobian transpose */ 2813b1b0ffSHong Zhang Mat MatDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 2913b1b0ffSHong Zhang Vec VecDeltaFwdSensipTemp; /* Working vector for holding one column of the sensitivity matrix */ 3013b1b0ffSHong Zhang Vec VecDeltaFwdSensipTemp2; /* Additional working vector for endpoint case */ 3113b1b0ffSHong Zhang Mat MatFwdSensip0; /* backup for roll-backs due to events */ 32715f1b00SHong Zhang Vec VecIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 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 170470880abSPatrick Sanan static PetscErrorCode TSTheta_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); 216470880abSPatrick Sanan ierr = TSTheta_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) { 278*c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 27905755b9cSHong Zhang } else { 280*c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);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) { 287*c9ad9de0SHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[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) { 311*c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);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) { 316*c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdu[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) { 327*c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 328b5e0532dSHong Zhang } 329b5e0532dSHong Zhang } 330b5e0532dSHong Zhang } 3313fd52205SHong Zhang 3323fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 333*c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(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.) { 339*c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(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) { 346*c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(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.) { 351*c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(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) { 362*c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);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) { 368*c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 36936eaed60SHong Zhang } 3702ca6e920SHong Zhang } 3713fd52205SHong Zhang if (ts->vecs_sensip) { 372*c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(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) { 378*c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(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; 43313b1b0ffSHong Zhang PetscInt 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; 44213b1b0ffSHong Zhang if (ts->mat_sensip) { 44313b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 4446f92202bSHong Zhang } 4456f92202bSHong Zhang if (ts->vecs_integral_sensip) { 4466f92202bSHong Zhang for (ncost=0;ncost<ts->numcost;ncost++) { 4476f92202bSHong Zhang ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 4486f92202bSHong Zhang } 4496f92202bSHong Zhang } 450715f1b00SHong Zhang PetscFunctionReturn(0); 451715f1b00SHong Zhang } 452715f1b00SHong Zhang 453715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 454715f1b00SHong Zhang { 455715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 45613b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 45713b1b0ffSHong Zhang Vec VecDeltaFwdSensipTemp = th->VecDeltaFwdSensipTemp,VecDeltaFwdSensipTemp2 = th->VecDeltaFwdSensipTemp2; 45813b1b0ffSHong Zhang PetscInt ncost,ntlm; 459715f1b00SHong Zhang KSP ksp; 460715f1b00SHong Zhang Mat J,Jp; 461715f1b00SHong Zhang PetscReal shift; 46213b1b0ffSHong Zhang PetscScalar *barr,*xarr; 463715f1b00SHong Zhang PetscErrorCode ierr; 464715f1b00SHong Zhang 465715f1b00SHong Zhang PetscFunctionBegin; 46613b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 46713b1b0ffSHong Zhang 4686f92202bSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 4696f92202bSHong Zhang if (ts->vecs_integral_sensip) { 4706f92202bSHong Zhang ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 4716f92202bSHong Zhang } 4726f92202bSHong Zhang } 4736f92202bSHong Zhang 474715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 475715f1b00SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 476715f1b00SHong Zhang 477715f1b00SHong Zhang /* Build RHS */ 478715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 479715f1b00SHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 480715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 48113b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 48213b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 483715f1b00SHong Zhang 484022c081aSHong Zhang /* Add the f_p forcing terms */ 485*c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 48613b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(1.-th->Theta)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 487*c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 48813b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 489715f1b00SHong Zhang } else { /* 1-stage method */ 490715f1b00SHong Zhang shift = 0.0; 491715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 49213b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 49313b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 49413b1b0ffSHong Zhang 495022c081aSHong Zhang /* Add the f_p forcing terms */ 496*c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 49713b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 498715f1b00SHong Zhang } 499715f1b00SHong Zhang 500715f1b00SHong Zhang /* Build LHS */ 501715f1b00SHong Zhang shift = 1/(th->Theta*th->time_step); 502715f1b00SHong Zhang if (th->endpoint) { 503715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 504715f1b00SHong Zhang } else { 505715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 506715f1b00SHong Zhang } 507715f1b00SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 508715f1b00SHong Zhang 509715f1b00SHong Zhang /* 510715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 511*c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 512715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 513715f1b00SHong Zhang */ 514715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 515715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 516*c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 517*c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 518715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 519*c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 52013b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 521715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 522715f1b00SHong Zhang } 523715f1b00SHong Zhang } 524715f1b00SHong Zhang } 525715f1b00SHong Zhang 526715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 527022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 52813b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 52913b1b0ffSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipTemp,barr);CHKERRQ(ierr); 530715f1b00SHong Zhang if (th->endpoint) { 53113b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 53213b1b0ffSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipTemp2,xarr);CHKERRQ(ierr); 53313b1b0ffSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipTemp,VecDeltaFwdSensipTemp2);CHKERRQ(ierr); 53413b1b0ffSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipTemp2);CHKERRQ(ierr); 53513b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 536715f1b00SHong Zhang } else { 53713b1b0ffSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipTemp,VecDeltaFwdSensipTemp);CHKERRQ(ierr); 538715f1b00SHong Zhang } 53913b1b0ffSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipTemp);CHKERRQ(ierr); 54013b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 541715f1b00SHong Zhang } 542715f1b00SHong Zhang 54313b1b0ffSHong Zhang /* 54413b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 545*c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 54613b1b0ffSHong Zhang */ 54713b1b0ffSHong Zhang if (ts->vecs_integral_sensip) { 54813b1b0ffSHong Zhang if (!th->endpoint) { 54913b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 550*c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 551*c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 55213b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 553*c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 55413b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 55513b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 55613b1b0ffSHong Zhang } 55713b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 55813b1b0ffSHong Zhang } else { 559*c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 560*c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 56113b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 562*c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 56313b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 56413b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 56513b1b0ffSHong Zhang } 56613b1b0ffSHong Zhang } 56713b1b0ffSHong Zhang } else { 56813b1b0ffSHong Zhang if (!th->endpoint) { 56913b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 570715f1b00SHong Zhang } 571715f1b00SHong Zhang } 5721566a47fSLisandro Dalcin PetscFunctionReturn(0); 5731566a47fSLisandro Dalcin } 5741566a47fSLisandro Dalcin 575316643e7SJed Brown /*------------------------------------------------------------*/ 576277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 577316643e7SJed Brown { 578316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 579316643e7SJed Brown PetscErrorCode ierr; 580316643e7SJed Brown 581316643e7SJed Brown PetscFunctionBegin; 5826bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 5836bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 5843b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 585eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 5861566a47fSLisandro Dalcin 5871566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 5881566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 5891566a47fSLisandro Dalcin 590b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 591715f1b00SHong Zhang if (ts->forward_solve) { 592715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 593715f1b00SHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 5946f92202bSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 595715f1b00SHong Zhang } 59613b1b0ffSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipTemp);CHKERRQ(ierr); 59713b1b0ffSHong Zhang if (th->endpoint) { 59813b1b0ffSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipTemp2);CHKERRQ(ierr); 59913b1b0ffSHong Zhang } 60013b1b0ffSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 60113b1b0ffSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 602715f1b00SHong Zhang } 603b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 604b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 605b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 606715f1b00SHong Zhang 607277b19d0SLisandro Dalcin PetscFunctionReturn(0); 608277b19d0SLisandro Dalcin } 609277b19d0SLisandro Dalcin 610277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 611277b19d0SLisandro Dalcin { 612277b19d0SLisandro Dalcin PetscErrorCode ierr; 613277b19d0SLisandro Dalcin 614277b19d0SLisandro Dalcin PetscFunctionBegin; 615277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 616b3a6b972SJed Brown if (ts->dm) { 617b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 618b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 619b3a6b972SJed Brown } 620277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 621bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 622bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 623bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 624bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 625316643e7SJed Brown PetscFunctionReturn(0); 626316643e7SJed Brown } 627316643e7SJed Brown 628316643e7SJed Brown /* 629316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 6302b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 631316643e7SJed Brown */ 6320f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 633316643e7SJed Brown { 634316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 635316643e7SJed Brown PetscErrorCode ierr; 6367445fe48SJed Brown Vec X0,Xdot; 6377445fe48SJed Brown DM dm,dmsave; 6381566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 639316643e7SJed Brown 640316643e7SJed Brown PetscFunctionBegin; 6417445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 6425a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 6437445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 644b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 6457445fe48SJed Brown 6467445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 6477445fe48SJed Brown dmsave = ts->dm; 6487445fe48SJed Brown ts->dm = dm; 6497445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 6507445fe48SJed Brown ts->dm = dmsave; 6510d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 652316643e7SJed Brown PetscFunctionReturn(0); 653316643e7SJed Brown } 654316643e7SJed Brown 655d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 656316643e7SJed Brown { 657316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 658316643e7SJed Brown PetscErrorCode ierr; 6597445fe48SJed Brown Vec Xdot; 6607445fe48SJed Brown DM dm,dmsave; 6611566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 662316643e7SJed Brown 663316643e7SJed Brown PetscFunctionBegin; 6647445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 665be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 6660298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 6677445fe48SJed Brown 6687445fe48SJed Brown dmsave = ts->dm; 6697445fe48SJed Brown ts->dm = dm; 670d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 6717445fe48SJed Brown ts->dm = dmsave; 6720298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 673316643e7SJed Brown PetscFunctionReturn(0); 674316643e7SJed Brown } 675316643e7SJed Brown 676715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 677715f1b00SHong Zhang { 678715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 679715f1b00SHong Zhang PetscErrorCode ierr; 680715f1b00SHong Zhang 681715f1b00SHong Zhang PetscFunctionBegin; 682022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 68313b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 68413b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 685715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 686715f1b00SHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 687715f1b00SHong Zhang } 6886f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 68913b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 69013b1b0ffSHong Zhang 6916f92202bSHong Zhang if (ts->vecs_integral_sensip) { 6926f92202bSHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 6936f92202bSHong Zhang } 69413b1b0ffSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipTemp);CHKERRQ(ierr); 69513b1b0ffSHong Zhang if (th->endpoint) { 69613b1b0ffSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipTemp2);CHKERRQ(ierr); 69713b1b0ffSHong Zhang } 698715f1b00SHong Zhang PetscFunctionReturn(0); 699715f1b00SHong Zhang } 700715f1b00SHong Zhang 701316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 702316643e7SJed Brown { 703316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 7042ffb9264SLisandro Dalcin PetscBool match; 705316643e7SJed Brown PetscErrorCode ierr; 706316643e7SJed Brown 707316643e7SJed Brown PetscFunctionBegin; 708b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 709b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 710b1cb36f3SHong Zhang } 71139d13666SHong Zhang if (!th->X) { 712316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 71339d13666SHong Zhang } 71439d13666SHong Zhang if (!th->Xdot) { 715316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 71639d13666SHong Zhang } 71739d13666SHong Zhang if (!th->X0) { 7183b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 71939d13666SHong Zhang } 7201566a47fSLisandro Dalcin if (th->endpoint) { 7211566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 7227445fe48SJed Brown } 7233b1890cdSShri Abhyankar 7241566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 7251566a47fSLisandro Dalcin 7261566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 7271566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 7281566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 7291566a47fSLisandro Dalcin 7301566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 7311566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 7322ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 7332ffb9264SLisandro Dalcin if (!match) { 7341566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 7351566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 7363b1890cdSShri Abhyankar } 7371566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 738b4dd3bf9SBarry Smith PetscFunctionReturn(0); 739b4dd3bf9SBarry Smith } 7400c7376e5SHong Zhang 741b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 742b4dd3bf9SBarry Smith 743b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 744b4dd3bf9SBarry Smith { 745b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 746b4dd3bf9SBarry Smith PetscErrorCode ierr; 747b4dd3bf9SBarry Smith 748b4dd3bf9SBarry Smith PetscFunctionBegin; 749b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 7502ca6e920SHong Zhang if(ts->vecs_sensip) { 751b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 7522ca6e920SHong Zhang } 753b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 754316643e7SJed Brown PetscFunctionReturn(0); 755316643e7SJed Brown } 756316643e7SJed Brown 7574416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 758316643e7SJed Brown { 759316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 760316643e7SJed Brown PetscErrorCode ierr; 761316643e7SJed Brown 762316643e7SJed Brown PetscFunctionBegin; 763e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 764316643e7SJed Brown { 7650298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 7660298fd71SBarry 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); 76703842d09SLisandro 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); 768316643e7SJed Brown } 769316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 770316643e7SJed Brown PetscFunctionReturn(0); 771316643e7SJed Brown } 772316643e7SJed Brown 773316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 774316643e7SJed Brown { 775316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 776ace3abfcSBarry Smith PetscBool iascii; 777316643e7SJed Brown PetscErrorCode ierr; 778316643e7SJed Brown 779316643e7SJed Brown PetscFunctionBegin; 780251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 781316643e7SJed Brown if (iascii) { 7827c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 783ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 784316643e7SJed Brown } 785316643e7SJed Brown PetscFunctionReturn(0); 786316643e7SJed Brown } 787316643e7SJed Brown 788560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 7890de4c49aSJed Brown { 7900de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 7910de4c49aSJed Brown 7920de4c49aSJed Brown PetscFunctionBegin; 7930de4c49aSJed Brown *theta = th->Theta; 7940de4c49aSJed Brown PetscFunctionReturn(0); 7950de4c49aSJed Brown } 7960de4c49aSJed Brown 797560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 7980de4c49aSJed Brown { 7990de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 8000de4c49aSJed Brown 8010de4c49aSJed Brown PetscFunctionBegin; 8027c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 8030de4c49aSJed Brown th->Theta = theta; 8041566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 8050de4c49aSJed Brown PetscFunctionReturn(0); 8060de4c49aSJed Brown } 807eb284becSJed Brown 808560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 80926f2ff8fSLisandro Dalcin { 81026f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 81126f2ff8fSLisandro Dalcin 81226f2ff8fSLisandro Dalcin PetscFunctionBegin; 81326f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 81426f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 81526f2ff8fSLisandro Dalcin } 81626f2ff8fSLisandro Dalcin 817560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 818eb284becSJed Brown { 819eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 820eb284becSJed Brown 821eb284becSJed Brown PetscFunctionBegin; 822eb284becSJed Brown th->endpoint = flg; 823eb284becSJed Brown PetscFunctionReturn(0); 824eb284becSJed Brown } 8250de4c49aSJed Brown 826f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 827f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 828f9c1d6abSBarry Smith { 829f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 830f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 8313fd8ae06SJed Brown const PetscReal one = 1.0; 832f9c1d6abSBarry Smith 833f9c1d6abSBarry Smith PetscFunctionBegin; 8343fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 835f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 836f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 837f9c1d6abSBarry Smith PetscFunctionReturn(0); 838f9c1d6abSBarry Smith } 839f9c1d6abSBarry Smith #endif 840f9c1d6abSBarry Smith 84142682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 84242682096SHong Zhang { 84342682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 84442682096SHong Zhang 84542682096SHong Zhang PetscFunctionBegin; 8461566a47fSLisandro Dalcin if (ns) *ns = 1; 8471566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 84842682096SHong Zhang PetscFunctionReturn(0); 84942682096SHong Zhang } 850f9c1d6abSBarry Smith 851316643e7SJed Brown /* ------------------------------------------------------------ */ 852316643e7SJed Brown /*MC 85396f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 854316643e7SJed Brown 855316643e7SJed Brown Level: beginner 856316643e7SJed Brown 8574eb428fdSBarry Smith Options Database: 858c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 859c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 86003842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 8614eb428fdSBarry Smith 862eb284becSJed Brown Notes: 8630c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 8640c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 8654eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 8664eb428fdSBarry Smith 867eb284becSJed Brown This method can be applied to DAE. 868eb284becSJed Brown 869eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 870eb284becSJed Brown 871eb284becSJed Brown .vb 872eb284becSJed Brown Theta | Theta 873eb284becSJed Brown ------------- 874eb284becSJed Brown | 1 875eb284becSJed Brown .ve 876eb284becSJed Brown 877eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 878eb284becSJed Brown 879eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 880eb284becSJed Brown 881eb284becSJed Brown .vb 882eb284becSJed Brown 0 | 0 0 883eb284becSJed Brown 1 | 1-Theta Theta 884eb284becSJed Brown ------------------- 885eb284becSJed Brown | 1-Theta Theta 886eb284becSJed Brown .ve 887eb284becSJed Brown 888eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 889eb284becSJed Brown 890eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 891eb284becSJed Brown 892eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 893eb284becSJed Brown 8944eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 895eb284becSJed Brown 896eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 897316643e7SJed Brown 898316643e7SJed Brown M*/ 8998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 900316643e7SJed Brown { 901316643e7SJed Brown TS_Theta *th; 902316643e7SJed Brown PetscErrorCode ierr; 903316643e7SJed Brown 904316643e7SJed Brown PetscFunctionBegin; 905277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 906316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 907316643e7SJed Brown ts->ops->view = TSView_Theta; 908316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 90942f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 910316643e7SJed Brown ts->ops->step = TSStep_Theta; 911cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 9121566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 91324655328SShri ts->ops->rollback = TSRollBack_Theta; 914316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 9150f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 9160f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 917f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 918f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 919f9c1d6abSBarry Smith #endif 92042682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 92142f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 922b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 923b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 9242ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 925715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 926715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 927316643e7SJed Brown 928efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 929efd4aadfSBarry Smith 930b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 931316643e7SJed Brown ts->data = (void*)th; 932316643e7SJed Brown 933715f1b00SHong Zhang th->VecsDeltaLam = NULL; 934715f1b00SHong Zhang th->VecsDeltaMu = NULL; 935715f1b00SHong Zhang th->VecsSensiTemp = NULL; 936715f1b00SHong Zhang 9376f700aefSJed Brown th->extrapolate = PETSC_FALSE; 938316643e7SJed Brown th->Theta = 0.5; 9391566a47fSLisandro Dalcin th->order = 2; 940bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 941bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 942bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 943bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 944316643e7SJed Brown PetscFunctionReturn(0); 945316643e7SJed Brown } 9460de4c49aSJed Brown 9470de4c49aSJed Brown /*@ 9480de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 9490de4c49aSJed Brown 9500de4c49aSJed Brown Not Collective 9510de4c49aSJed Brown 9520de4c49aSJed Brown Input Parameter: 9530de4c49aSJed Brown . ts - timestepping context 9540de4c49aSJed Brown 9550de4c49aSJed Brown Output Parameter: 9560de4c49aSJed Brown . theta - stage abscissa 9570de4c49aSJed Brown 9580de4c49aSJed Brown Note: 9590de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 9600de4c49aSJed Brown 9610de4c49aSJed Brown Level: Advanced 9620de4c49aSJed Brown 9630de4c49aSJed Brown .seealso: TSThetaSetTheta() 9640de4c49aSJed Brown @*/ 9657087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 9660de4c49aSJed Brown { 9674ac538c5SBarry Smith PetscErrorCode ierr; 9680de4c49aSJed Brown 9690de4c49aSJed Brown PetscFunctionBegin; 970afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 9710de4c49aSJed Brown PetscValidPointer(theta,2); 9724ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 9730de4c49aSJed Brown PetscFunctionReturn(0); 9740de4c49aSJed Brown } 9750de4c49aSJed Brown 9760de4c49aSJed Brown /*@ 9770de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 9780de4c49aSJed Brown 9790de4c49aSJed Brown Not Collective 9800de4c49aSJed Brown 9810de4c49aSJed Brown Input Parameter: 9820de4c49aSJed Brown + ts - timestepping context 9830de4c49aSJed Brown - theta - stage abscissa 9840de4c49aSJed Brown 9850de4c49aSJed Brown Options Database: 9860de4c49aSJed Brown . -ts_theta_theta <theta> 9870de4c49aSJed Brown 9880de4c49aSJed Brown Level: Intermediate 9890de4c49aSJed Brown 9900de4c49aSJed Brown .seealso: TSThetaGetTheta() 9910de4c49aSJed Brown @*/ 9927087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 9930de4c49aSJed Brown { 9944ac538c5SBarry Smith PetscErrorCode ierr; 9950de4c49aSJed Brown 9960de4c49aSJed Brown PetscFunctionBegin; 997afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 9984ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 9990de4c49aSJed Brown PetscFunctionReturn(0); 10000de4c49aSJed Brown } 1001f33bbcb6SJed Brown 100226f2ff8fSLisandro Dalcin /*@ 100326f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 100426f2ff8fSLisandro Dalcin 100526f2ff8fSLisandro Dalcin Not Collective 100626f2ff8fSLisandro Dalcin 100726f2ff8fSLisandro Dalcin Input Parameter: 100826f2ff8fSLisandro Dalcin . ts - timestepping context 100926f2ff8fSLisandro Dalcin 101026f2ff8fSLisandro Dalcin Output Parameter: 101126f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 101226f2ff8fSLisandro Dalcin 101326f2ff8fSLisandro Dalcin Level: Advanced 101426f2ff8fSLisandro Dalcin 101526f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 101626f2ff8fSLisandro Dalcin @*/ 101726f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 101826f2ff8fSLisandro Dalcin { 101926f2ff8fSLisandro Dalcin PetscErrorCode ierr; 102026f2ff8fSLisandro Dalcin 102126f2ff8fSLisandro Dalcin PetscFunctionBegin; 102226f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 102326f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1024163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 102526f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 102626f2ff8fSLisandro Dalcin } 102726f2ff8fSLisandro Dalcin 1028eb284becSJed Brown /*@ 1029eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1030eb284becSJed Brown 1031eb284becSJed Brown Not Collective 1032eb284becSJed Brown 1033eb284becSJed Brown Input Parameter: 1034eb284becSJed Brown + ts - timestepping context 1035eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1036eb284becSJed Brown 1037eb284becSJed Brown Options Database: 1038eb284becSJed Brown . -ts_theta_endpoint <flg> 1039eb284becSJed Brown 1040eb284becSJed Brown Level: Intermediate 1041eb284becSJed Brown 1042eb284becSJed Brown .seealso: TSTHETA, TSCN 1043eb284becSJed Brown @*/ 1044eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1045eb284becSJed Brown { 1046eb284becSJed Brown PetscErrorCode ierr; 1047eb284becSJed Brown 1048eb284becSJed Brown PetscFunctionBegin; 1049eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1050eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1051eb284becSJed Brown PetscFunctionReturn(0); 1052eb284becSJed Brown } 1053eb284becSJed Brown 1054f33bbcb6SJed Brown /* 1055f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1056f33bbcb6SJed Brown * The creation functions for these specializations are below. 1057f33bbcb6SJed Brown */ 1058f33bbcb6SJed Brown 10591566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 10601566a47fSLisandro Dalcin { 10611566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 10621566a47fSLisandro Dalcin PetscErrorCode ierr; 10631566a47fSLisandro Dalcin 10641566a47fSLisandro Dalcin PetscFunctionBegin; 10651566a47fSLisandro 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"); 10661566a47fSLisandro 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"); 10671566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 10681566a47fSLisandro Dalcin PetscFunctionReturn(0); 10691566a47fSLisandro Dalcin } 10701566a47fSLisandro Dalcin 1071f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1072f33bbcb6SJed Brown { 1073f33bbcb6SJed Brown PetscFunctionBegin; 1074f33bbcb6SJed Brown PetscFunctionReturn(0); 1075f33bbcb6SJed Brown } 1076f33bbcb6SJed Brown 1077f33bbcb6SJed Brown /*MC 1078f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1079f33bbcb6SJed Brown 1080f33bbcb6SJed Brown Level: beginner 1081f33bbcb6SJed Brown 10824eb428fdSBarry Smith Notes: 1083c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 10844eb428fdSBarry Smith 10851566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 10864eb428fdSBarry Smith 1087f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1088f33bbcb6SJed Brown 1089f33bbcb6SJed Brown M*/ 10908cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1091f33bbcb6SJed Brown { 1092f33bbcb6SJed Brown PetscErrorCode ierr; 1093f33bbcb6SJed Brown 1094f33bbcb6SJed Brown PetscFunctionBegin; 1095f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1096f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 10971566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 10980c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1099f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1100f33bbcb6SJed Brown PetscFunctionReturn(0); 1101f33bbcb6SJed Brown } 1102f33bbcb6SJed Brown 11031566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 11041566a47fSLisandro Dalcin { 11051566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 11061566a47fSLisandro Dalcin PetscErrorCode ierr; 11071566a47fSLisandro Dalcin 11081566a47fSLisandro Dalcin PetscFunctionBegin; 11091566a47fSLisandro 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"); 11101566a47fSLisandro 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"); 11111566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 11121566a47fSLisandro Dalcin PetscFunctionReturn(0); 11131566a47fSLisandro Dalcin } 11141566a47fSLisandro Dalcin 1115f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1116f33bbcb6SJed Brown { 1117f33bbcb6SJed Brown PetscFunctionBegin; 1118f33bbcb6SJed Brown PetscFunctionReturn(0); 1119f33bbcb6SJed Brown } 1120f33bbcb6SJed Brown 1121f33bbcb6SJed Brown /*MC 1122f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1123f33bbcb6SJed Brown 1124f33bbcb6SJed Brown Level: beginner 1125f33bbcb6SJed Brown 1126f33bbcb6SJed Brown Notes: 11277cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 11287cf5af47SJed Brown 11297cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1130f33bbcb6SJed Brown 1131f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1132f33bbcb6SJed Brown 1133f33bbcb6SJed Brown M*/ 11348cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1135f33bbcb6SJed Brown { 1136f33bbcb6SJed Brown PetscErrorCode ierr; 1137f33bbcb6SJed Brown 1138f33bbcb6SJed Brown PetscFunctionBegin; 1139f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1140f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1141eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 11420c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1143f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1144f33bbcb6SJed Brown PetscFunctionReturn(0); 1145f33bbcb6SJed Brown } 1146