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 */ 29b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */ 3013b1b0ffSHong Zhang Mat MatFwdSensip0; /* backup for roll-backs due to events */ 31715f1b00SHong Zhang Vec VecIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 326f92202bSHong Zhang Vec *VecsIntegralSensip0; /* backup for roll-backs due to events */ 33b5b4017aSHong Zhang Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */ 34b5b4017aSHong Zhang Vec *VecsDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */ 35b5b4017aSHong Zhang Vec *VecsSensi2Temp; /* Working vectors that holds the residual for the second-order adjoint */ 36b5b4017aSHong Zhang Vec *VecsAffine; /* Working vectors to store residuals */ 37715f1b00SHong Zhang /* context for error estimation */ 381566a47fSLisandro Dalcin Vec vec_sol_prev; 391566a47fSLisandro Dalcin Vec vec_lte_work; 40316643e7SJed Brown } TS_Theta; 41316643e7SJed Brown 427445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 437445fe48SJed Brown { 447445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 457445fe48SJed Brown PetscErrorCode ierr; 467445fe48SJed Brown 477445fe48SJed Brown PetscFunctionBegin; 487445fe48SJed Brown if (X0) { 497445fe48SJed Brown if (dm && dm != ts->dm) { 500d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 517445fe48SJed Brown } else *X0 = ts->vec_sol; 527445fe48SJed Brown } 537445fe48SJed Brown if (Xdot) { 547445fe48SJed Brown if (dm && dm != ts->dm) { 550d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 567445fe48SJed Brown } else *Xdot = th->Xdot; 577445fe48SJed Brown } 587445fe48SJed Brown PetscFunctionReturn(0); 597445fe48SJed Brown } 607445fe48SJed Brown 610d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 620d0b770aSPeter Brune { 630d0b770aSPeter Brune PetscErrorCode ierr; 640d0b770aSPeter Brune 650d0b770aSPeter Brune PetscFunctionBegin; 660d0b770aSPeter Brune if (X0) { 670d0b770aSPeter Brune if (dm && dm != ts->dm) { 680d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 690d0b770aSPeter Brune } 700d0b770aSPeter Brune } 710d0b770aSPeter Brune if (Xdot) { 720d0b770aSPeter Brune if (dm && dm != ts->dm) { 730d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 740d0b770aSPeter Brune } 750d0b770aSPeter Brune } 760d0b770aSPeter Brune PetscFunctionReturn(0); 770d0b770aSPeter Brune } 780d0b770aSPeter Brune 797445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 807445fe48SJed Brown { 817445fe48SJed Brown PetscFunctionBegin; 827445fe48SJed Brown PetscFunctionReturn(0); 837445fe48SJed Brown } 847445fe48SJed Brown 857445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 867445fe48SJed Brown { 877445fe48SJed Brown TS ts = (TS)ctx; 887445fe48SJed Brown PetscErrorCode ierr; 897445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 907445fe48SJed Brown 917445fe48SJed Brown PetscFunctionBegin; 927445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 937445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 947445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 957445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 967445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 977445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 980d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 990d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 1007445fe48SJed Brown PetscFunctionReturn(0); 1017445fe48SJed Brown } 1027445fe48SJed Brown 103258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 104258e1594SPeter Brune { 105258e1594SPeter Brune PetscFunctionBegin; 106258e1594SPeter Brune PetscFunctionReturn(0); 107258e1594SPeter Brune } 108258e1594SPeter Brune 109258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 110258e1594SPeter Brune { 111258e1594SPeter Brune TS ts = (TS)ctx; 112258e1594SPeter Brune PetscErrorCode ierr; 113258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 114258e1594SPeter Brune 115258e1594SPeter Brune PetscFunctionBegin; 116258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 117258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 118258e1594SPeter Brune 119258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 120258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121258e1594SPeter Brune 122258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 123258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 124258e1594SPeter Brune 125258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 126258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 127258e1594SPeter Brune PetscFunctionReturn(0); 128258e1594SPeter Brune } 129258e1594SPeter Brune 130022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 131022c081aSHong Zhang { 132022c081aSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 133022c081aSHong Zhang PetscErrorCode ierr; 134022c081aSHong Zhang 135022c081aSHong Zhang PetscFunctionBegin; 136022c081aSHong Zhang if (th->endpoint) { 137022c081aSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 138022c081aSHong Zhang if (th->Theta!=1.0) { 139022c081aSHong Zhang ierr = TSComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 140022c081aSHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 141022c081aSHong Zhang } 142022c081aSHong Zhang ierr = TSComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 143022c081aSHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 144022c081aSHong Zhang } else { 145022c081aSHong Zhang ierr = TSComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 146022c081aSHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 147022c081aSHong Zhang } 148022c081aSHong Zhang PetscFunctionReturn(0); 149022c081aSHong Zhang } 150022c081aSHong Zhang 151b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 152b1cb36f3SHong Zhang { 153b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 154b1cb36f3SHong Zhang PetscErrorCode ierr; 155b1cb36f3SHong Zhang 156b1cb36f3SHong Zhang PetscFunctionBegin; 157b1cb36f3SHong Zhang /* backup cost integral */ 158b1cb36f3SHong Zhang ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr); 159022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 160b1cb36f3SHong Zhang PetscFunctionReturn(0); 161b1cb36f3SHong Zhang } 162b1cb36f3SHong Zhang 163b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 164b1cb36f3SHong Zhang { 165b1cb36f3SHong Zhang PetscErrorCode ierr; 166b1cb36f3SHong Zhang 167b1cb36f3SHong Zhang PetscFunctionBegin; 168022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 16924655328SShri PetscFunctionReturn(0); 17024655328SShri } 17124655328SShri 172470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x) 1731566a47fSLisandro Dalcin { 1741566a47fSLisandro Dalcin PetscInt nits,lits; 1751566a47fSLisandro Dalcin PetscErrorCode ierr; 1761566a47fSLisandro Dalcin 1771566a47fSLisandro Dalcin PetscFunctionBegin; 1781566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1791566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1801566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1811566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1821566a47fSLisandro Dalcin PetscFunctionReturn(0); 1831566a47fSLisandro Dalcin } 1841566a47fSLisandro Dalcin 185193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 186316643e7SJed Brown { 187316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1881566a47fSLisandro Dalcin PetscInt rejections = 0; 1894957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 1901566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 191051f2191SLisandro Dalcin PetscErrorCode ierr; 192316643e7SJed Brown 193316643e7SJed Brown PetscFunctionBegin; 1941566a47fSLisandro Dalcin if (!ts->steprollback) { 1952ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 1963b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 1971566a47fSLisandro Dalcin } 1981566a47fSLisandro Dalcin 1991566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 2001566a47fSLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 201715f1b00SHong Zhang 2021566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 2031566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 204316643e7SJed Brown 205be5899b3SLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 206fecfb714SLisandro Dalcin if (th->extrapolate && !ts->steprestart) { 207be5899b3SLisandro Dalcin ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 208be5899b3SLisandro Dalcin } 209eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 210eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2111566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2121566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2131566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2141566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2151566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 216eb284becSJed Brown } 217be5899b3SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 218470880abSPatrick Sanan ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 219be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 2201566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 221fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 222051f2191SLisandro Dalcin 223051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2241566a47fSLisandro Dalcin if (th->endpoint) { 2251566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2261566a47fSLisandro Dalcin } else { 227cb3a5882SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 2281566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2291566a47fSLisandro Dalcin } 2301566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2311566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2321566a47fSLisandro Dalcin if (!accept) { 2331566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 234be5899b3SLisandro Dalcin ts->time_step = next_time_step; 235be5899b3SLisandro Dalcin goto reject_step; 236051f2191SLisandro Dalcin } 2373b1890cdSShri Abhyankar 238715f1b00SHong Zhang if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 239b1cb36f3SHong Zhang th->ptime = ts->ptime; 240b1cb36f3SHong Zhang th->time_step = ts->time_step; 24117145e75SHong Zhang } 2421566a47fSLisandro Dalcin 2432b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 244cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 245051f2191SLisandro Dalcin break; 246051f2191SLisandro Dalcin 247051f2191SLisandro Dalcin reject_step: 248fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 2491566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 250051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2511566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2523b1890cdSShri Abhyankar } 2533b1890cdSShri Abhyankar } 254316643e7SJed Brown PetscFunctionReturn(0); 255316643e7SJed Brown } 256316643e7SJed Brown 25742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 2582ca6e920SHong Zhang { 2592ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 260b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 261b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 2622ca6e920SHong Zhang PetscInt nadj; 2632ca6e920SHong Zhang Mat J,Jp; 2642ca6e920SHong Zhang KSP ksp; 2652ca6e920SHong Zhang PetscReal shift; 266b5b4017aSHong Zhang PetscScalar *xarr; 267b5b4017aSHong Zhang PetscErrorCode ierr; 2682ca6e920SHong Zhang 2692ca6e920SHong Zhang PetscFunctionBegin; 2702ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 271302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 2722ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 273b5e0532dSHong Zhang 274b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 275022c081aSHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 276b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 277022c081aSHong Zhang th->time_step = -ts->time_step; 2782ca6e920SHong Zhang 279b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 2802c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 28105755b9cSHong Zhang if (th->endpoint) { 282c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 28305755b9cSHong Zhang } else { 284c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 28505755b9cSHong Zhang } 28605755b9cSHong Zhang } 287abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 288b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 289022c081aSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr); 2902c39e106SBarry Smith if (ts->vec_costintegral) { 291c9ad9de0SHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 29236eaed60SHong Zhang } 2932ca6e920SHong Zhang } 2943c54f92cSHong Zhang 295b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 296022c081aSHong Zhang shift = 1./(th->Theta*th->time_step); 2973c54f92cSHong Zhang if (th->endpoint) { 298022c081aSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 2993c54f92cSHong Zhang } else { 3003c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3013c54f92cSHong Zhang } 3022ca6e920SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 3032ca6e920SHong Zhang 304b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 305abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 306b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 3072ca6e920SHong Zhang } 3083c54f92cSHong Zhang 309b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 310b5b4017aSHong Zhang if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 311b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 312b5b4017aSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 313b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 314b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 315b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 316b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 317b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 318b5b4017aSHong Zhang if (ts->vecs_fup) { 319b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 320b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 321b5b4017aSHong Zhang } 322b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 323b5b4017aSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 324b5b4017aSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],1./shift);CHKERRQ(ierr); 325b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 326b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 327b5b4017aSHong Zhang if (ts->vecs_fup) { 328b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 329b5b4017aSHong Zhang } 330b5b4017aSHong Zhang if (ts->vec_costintegral) { 331b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 332b5b4017aSHong Zhang } 333b5b4017aSHong Zhang } 334b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 335b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 336b5b4017aSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 337b5b4017aSHong Zhang } 338b5b4017aSHong Zhang } 339b5b4017aSHong Zhang 34036eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 341b5e0532dSHong Zhang if(th->endpoint) { /* two-stage case */ 342b5b4017aSHong Zhang if (th->Theta != 1.) { /* general case */ 343022c081aSHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 344b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 345b5b4017aSHong Zhang if (ts->vec_costintegral) { /* R_U at t_n */ 346c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 3478263dbbfSHong Zhang } 348abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 349b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 3503c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 351b5b4017aSHong Zhang if (ts->vec_costintegral) { 352b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1./shift,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 353b5b4017aSHong Zhang } 354b5b4017aSHong Zhang } 355b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* second-order */ 356b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 357b5b4017aSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 358b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 359b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 360b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 361b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 362b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 363b5b4017aSHong Zhang if (ts->vecs_fup) { 364b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 365b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 366b5b4017aSHong Zhang } 367b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 368b5b4017aSHong Zhang /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) R_U */ 369b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 370b5b4017aSHong Zhang ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr); 371b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 372b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 373b5b4017aSHong Zhang if (ts->vecs_fup) { 374b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fup[nadj]);CHKERRQ(ierr); 375b5b4017aSHong Zhang } 376b5b4017aSHong Zhang if (ts->vec_costintegral) { 377b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 378b5b4017aSHong Zhang } 379b5b4017aSHong Zhang } 3802ca6e920SHong Zhang } 381b5e0532dSHong Zhang } else { /* backward Euler */ 382b5e0532dSHong Zhang shift = 0.0; 383b5b4017aSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_u */ 384b5e0532dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 385b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 386022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 387b5b4017aSHong Zhang if (ts->vec_costintegral) { /* wrong? */ 388c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 389b5e0532dSHong Zhang } 390b5e0532dSHong Zhang } 391b5b4017aSHong Zhang if (ts->vecs_sensi2) { 392b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 393b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 394b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-th->time_step,VecsSensi2Temp[nadj]);CHKERRQ(ierr); 395b5b4017aSHong Zhang } 396b5b4017aSHong Zhang } 397b5e0532dSHong Zhang } 3983fd52205SHong Zhang 3993fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 400b5b4017aSHong Zhang /* U_{n+1} */ 401c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 402b5b4017aSHong Zhang if (ts->vec_costintegral) { 403b5b4017aSHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 404b5b4017aSHong Zhang } 405abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 406b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 407022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 4083fd52205SHong Zhang } 409b5b4017aSHong Zhang if (ts->vecs_sensip2) { /* second-order */ 410b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 411b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 412b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 413b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 414b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 415b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 416b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 417b5b4017aSHong Zhang if (ts->vecs_fpu) { 418b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 419b5b4017aSHong Zhang } 420b5b4017aSHong Zhang if (ts->vecs_fpp) { 421b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 422b5b4017aSHong Zhang } 423b5b4017aSHong Zhang if (ts->vec_costintegral) { 424b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 425b5b4017aSHong Zhang } 426b5b4017aSHong Zhang } 427b5b4017aSHong Zhang } 428b5b4017aSHong Zhang 429b5b4017aSHong Zhang /* U_s */ 430b5e0532dSHong Zhang if (th->Theta!=1.) { 431c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 432b5b4017aSHong Zhang if (ts->vec_costintegral) { 433b5b4017aSHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 434b5b4017aSHong Zhang } 435abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 436b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 437022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 438b5b4017aSHong Zhang if (ts->vecs_sensip2) { /* second-order */ 439b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 440b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 441b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 442b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 443b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 444b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 445b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 446b5b4017aSHong Zhang if (ts->vecs_fpu) { 447b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 448b5e0532dSHong Zhang } 449b5b4017aSHong Zhang if (ts->vecs_fpp) { 450b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 4513fd52205SHong Zhang } 4522c39e106SBarry Smith if (ts->vec_costintegral) { 453b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 45436eaed60SHong Zhang } 455b5b4017aSHong Zhang } 45636eaed60SHong Zhang } 45736eaed60SHong Zhang } 4583fd52205SHong Zhang } 459b5e0532dSHong Zhang } 4603fd52205SHong Zhang } else { /* one-stage case */ 4613c54f92cSHong Zhang shift = 0.0; 462a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 4632c39e106SBarry Smith if (ts->vec_costintegral) { 464c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 4653c54f92cSHong Zhang } 466abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 467b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 468022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 4692c39e106SBarry Smith if (ts->vec_costintegral) { 470c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 47136eaed60SHong Zhang } 4722ca6e920SHong Zhang } 4733fd52205SHong Zhang if (ts->vecs_sensip) { 474c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 475abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 476b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 477022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 4783fd52205SHong Zhang } 4792c39e106SBarry Smith if (ts->vec_costintegral) { 480c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 481abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 482022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 48336eaed60SHong Zhang } 48436eaed60SHong Zhang } 4853fd52205SHong Zhang } 4862ca6e920SHong Zhang } 4872ca6e920SHong Zhang 4882ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 4892ca6e920SHong Zhang PetscFunctionReturn(0); 4902ca6e920SHong Zhang } 4912ca6e920SHong Zhang 492cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 493cd652676SJed Brown { 494cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 495be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 496cd652676SJed Brown PetscErrorCode ierr; 497cd652676SJed Brown 498cd652676SJed Brown PetscFunctionBegin; 499a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 500be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 501be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 502cd652676SJed Brown PetscFunctionReturn(0); 503cd652676SJed Brown } 504cd652676SJed Brown 5051566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 5061566a47fSLisandro Dalcin { 5071566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 5081566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 5091566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 5107453f775SEmil Constantinescu PetscReal wltea,wlter; 5111566a47fSLisandro Dalcin PetscErrorCode ierr; 5121566a47fSLisandro Dalcin 5131566a47fSLisandro Dalcin PetscFunctionBegin; 5142ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 5151566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 516fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 5171566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 518fecfb714SLisandro Dalcin { 519be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 520be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 5211566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 5221566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 5231566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 5241566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 5251566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 5267453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 5271566a47fSLisandro Dalcin } 5281566a47fSLisandro Dalcin if (order) *order = 2; 5291566a47fSLisandro Dalcin PetscFunctionReturn(0); 5301566a47fSLisandro Dalcin } 5311566a47fSLisandro Dalcin 5321566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 5331566a47fSLisandro Dalcin { 5341566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 53513b1b0ffSHong Zhang PetscInt ncost; 5361566a47fSLisandro Dalcin PetscErrorCode ierr; 5371566a47fSLisandro Dalcin 5381566a47fSLisandro Dalcin PetscFunctionBegin; 5391566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 5401566a47fSLisandro Dalcin if (ts->vec_costintegral && ts->costintegralfwd) { 5411566a47fSLisandro Dalcin ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 5421566a47fSLisandro Dalcin } 543715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 54413b1b0ffSHong Zhang if (ts->mat_sensip) { 54513b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 5466f92202bSHong Zhang } 5476f92202bSHong Zhang if (ts->vecs_integral_sensip) { 5486f92202bSHong Zhang for (ncost=0;ncost<ts->numcost;ncost++) { 5496f92202bSHong Zhang ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 5506f92202bSHong Zhang } 5516f92202bSHong Zhang } 552715f1b00SHong Zhang PetscFunctionReturn(0); 553715f1b00SHong Zhang } 554715f1b00SHong Zhang 555715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 556715f1b00SHong Zhang { 557715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 55813b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 559b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 56013b1b0ffSHong Zhang PetscInt ncost,ntlm; 561715f1b00SHong Zhang KSP ksp; 562715f1b00SHong Zhang Mat J,Jp; 563715f1b00SHong Zhang PetscReal shift; 56413b1b0ffSHong Zhang PetscScalar *barr,*xarr; 565715f1b00SHong Zhang PetscErrorCode ierr; 566715f1b00SHong Zhang 567715f1b00SHong Zhang PetscFunctionBegin; 56813b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 56913b1b0ffSHong Zhang 5706f92202bSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 5716f92202bSHong Zhang if (ts->vecs_integral_sensip) { 5726f92202bSHong Zhang ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 5736f92202bSHong Zhang } 5746f92202bSHong Zhang } 5756f92202bSHong Zhang 576715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 577715f1b00SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 578715f1b00SHong Zhang 579715f1b00SHong Zhang /* Build RHS */ 580715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 581715f1b00SHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 582715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 58313b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 58413b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 585715f1b00SHong Zhang 586022c081aSHong Zhang /* Add the f_p forcing terms */ 587*0e11c21fSHong Zhang if (ts->Jacp) { 588c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 58913b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(1.-th->Theta)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 590c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 59113b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 592*0e11c21fSHong Zhang } 593715f1b00SHong Zhang } else { /* 1-stage method */ 594715f1b00SHong Zhang shift = 0.0; 595715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 59613b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 59713b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 59813b1b0ffSHong Zhang 599022c081aSHong Zhang /* Add the f_p forcing terms */ 600*0e11c21fSHong Zhang if (ts->Jacp) { 601c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 60213b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 603715f1b00SHong Zhang } 604*0e11c21fSHong Zhang } 605715f1b00SHong Zhang 606715f1b00SHong Zhang /* Build LHS */ 607715f1b00SHong Zhang shift = 1/(th->Theta*th->time_step); 608715f1b00SHong Zhang if (th->endpoint) { 609715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 610715f1b00SHong Zhang } else { 611715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 612715f1b00SHong Zhang } 613715f1b00SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 614715f1b00SHong Zhang 615715f1b00SHong Zhang /* 616715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 617c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 618715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 619715f1b00SHong Zhang */ 620715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 621715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 622c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 623*0e11c21fSHong Zhang if (ts->vecs_drdp) { 624c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 625*0e11c21fSHong Zhang } 626715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 627c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 628*0e11c21fSHong Zhang if (ts->vecs_drdp) { 62913b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 630*0e11c21fSHong Zhang } 631715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 632715f1b00SHong Zhang } 633715f1b00SHong Zhang } 634715f1b00SHong Zhang } 635715f1b00SHong Zhang 636715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 637022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 638b5b4017aSHong Zhang KSPConvergedReason kspreason; 63913b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 640b5b4017aSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 641715f1b00SHong Zhang if (th->endpoint) { 64213b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 643b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 644b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 645b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 64613b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 647715f1b00SHong Zhang } else { 648b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 649715f1b00SHong Zhang } 650b5b4017aSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 651b5b4017aSHong Zhang if (kspreason < 0) { 652b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 653b5b4017aSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 654b5b4017aSHong Zhang } 655b5b4017aSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 65613b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 657715f1b00SHong Zhang } 658715f1b00SHong Zhang 659b5b4017aSHong Zhang 66013b1b0ffSHong Zhang /* 66113b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 662c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 66313b1b0ffSHong Zhang */ 66413b1b0ffSHong Zhang if (ts->vecs_integral_sensip) { 66513b1b0ffSHong Zhang if (!th->endpoint) { 66613b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 667c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 668*0e11c21fSHong Zhang if (ts->vecs_drdp) { 669c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 670*0e11c21fSHong Zhang } 67113b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 672c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 673*0e11c21fSHong Zhang if (ts->vecs_drdp) { 67413b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 675*0e11c21fSHong Zhang } 67613b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 67713b1b0ffSHong Zhang } 67813b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 67913b1b0ffSHong Zhang } else { 680c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 681*0e11c21fSHong Zhang if (ts->vecs_drdp) { 682c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 683*0e11c21fSHong Zhang } 68413b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 685c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 686*0e11c21fSHong Zhang if (ts->vecs_drdp) { 68713b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 688*0e11c21fSHong Zhang } 68913b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 69013b1b0ffSHong Zhang } 69113b1b0ffSHong Zhang } 69213b1b0ffSHong Zhang } else { 69313b1b0ffSHong Zhang if (!th->endpoint) { 69413b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 695715f1b00SHong Zhang } 696715f1b00SHong Zhang } 6971566a47fSLisandro Dalcin PetscFunctionReturn(0); 6981566a47fSLisandro Dalcin } 6991566a47fSLisandro Dalcin 700316643e7SJed Brown /*------------------------------------------------------------*/ 701277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 702316643e7SJed Brown { 703316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 704316643e7SJed Brown PetscErrorCode ierr; 705316643e7SJed Brown 706316643e7SJed Brown PetscFunctionBegin; 7076bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 7086bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 7093b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 710eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 7111566a47fSLisandro Dalcin 7121566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 7131566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 7141566a47fSLisandro Dalcin 715b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 716715f1b00SHong Zhang if (ts->forward_solve) { 717715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 718715f1b00SHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 7196f92202bSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 720715f1b00SHong Zhang } 721b5b4017aSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 72213b1b0ffSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 72313b1b0ffSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 724715f1b00SHong Zhang } 725b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 726b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 727b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 728b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 729b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 730b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 731715f1b00SHong Zhang 732277b19d0SLisandro Dalcin PetscFunctionReturn(0); 733277b19d0SLisandro Dalcin } 734277b19d0SLisandro Dalcin 735277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 736277b19d0SLisandro Dalcin { 737277b19d0SLisandro Dalcin PetscErrorCode ierr; 738277b19d0SLisandro Dalcin 739277b19d0SLisandro Dalcin PetscFunctionBegin; 740277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 741b3a6b972SJed Brown if (ts->dm) { 742b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 743b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 744b3a6b972SJed Brown } 745277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 746bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 747bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 748bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 749bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 750316643e7SJed Brown PetscFunctionReturn(0); 751316643e7SJed Brown } 752316643e7SJed Brown 753316643e7SJed Brown /* 754316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 7552b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 756316643e7SJed Brown */ 7570f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 758316643e7SJed Brown { 759316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 760316643e7SJed Brown PetscErrorCode ierr; 7617445fe48SJed Brown Vec X0,Xdot; 7627445fe48SJed Brown DM dm,dmsave; 7631566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 764316643e7SJed Brown 765316643e7SJed Brown PetscFunctionBegin; 7667445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 7675a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 7687445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 769b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 7707445fe48SJed Brown 7717445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 7727445fe48SJed Brown dmsave = ts->dm; 7737445fe48SJed Brown ts->dm = dm; 7747445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 7757445fe48SJed Brown ts->dm = dmsave; 7760d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 777316643e7SJed Brown PetscFunctionReturn(0); 778316643e7SJed Brown } 779316643e7SJed Brown 780d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 781316643e7SJed Brown { 782316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 783316643e7SJed Brown PetscErrorCode ierr; 7847445fe48SJed Brown Vec Xdot; 7857445fe48SJed Brown DM dm,dmsave; 7861566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 787316643e7SJed Brown 788316643e7SJed Brown PetscFunctionBegin; 7897445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 790be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 7910298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 7927445fe48SJed Brown 7937445fe48SJed Brown dmsave = ts->dm; 7947445fe48SJed Brown ts->dm = dm; 795d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 7967445fe48SJed Brown ts->dm = dmsave; 7970298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 798316643e7SJed Brown PetscFunctionReturn(0); 799316643e7SJed Brown } 800316643e7SJed Brown 801715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 802715f1b00SHong Zhang { 803715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 804715f1b00SHong Zhang PetscErrorCode ierr; 805715f1b00SHong Zhang 806715f1b00SHong Zhang PetscFunctionBegin; 807022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 80813b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 80913b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 810715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 811715f1b00SHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 812715f1b00SHong Zhang } 8136f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 81413b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 81513b1b0ffSHong Zhang 8166f92202bSHong Zhang if (ts->vecs_integral_sensip) { 8176f92202bSHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 8186f92202bSHong Zhang } 819b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 820b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 821715f1b00SHong Zhang PetscFunctionReturn(0); 822715f1b00SHong Zhang } 823715f1b00SHong Zhang 824316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 825316643e7SJed Brown { 826316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 8272ffb9264SLisandro Dalcin PetscBool match; 828316643e7SJed Brown PetscErrorCode ierr; 829316643e7SJed Brown 830316643e7SJed Brown PetscFunctionBegin; 831b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 832b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 833b1cb36f3SHong Zhang } 83439d13666SHong Zhang if (!th->X) { 835316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 83639d13666SHong Zhang } 83739d13666SHong Zhang if (!th->Xdot) { 838316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 83939d13666SHong Zhang } 84039d13666SHong Zhang if (!th->X0) { 8413b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 84239d13666SHong Zhang } 8431566a47fSLisandro Dalcin if (th->endpoint) { 8441566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 8457445fe48SJed Brown } 8463b1890cdSShri Abhyankar 8471566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 8481566a47fSLisandro Dalcin 8491566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 8501566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 8511566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 8521566a47fSLisandro Dalcin 8531566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 8541566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 8552ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 8562ffb9264SLisandro Dalcin if (!match) { 8571566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 8581566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 8593b1890cdSShri Abhyankar } 8601566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 861b4dd3bf9SBarry Smith PetscFunctionReturn(0); 862b4dd3bf9SBarry Smith } 8630c7376e5SHong Zhang 864b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 865b4dd3bf9SBarry Smith 866b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 867b4dd3bf9SBarry Smith { 868b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 869b4dd3bf9SBarry Smith PetscErrorCode ierr; 870b4dd3bf9SBarry Smith 871b4dd3bf9SBarry Smith PetscFunctionBegin; 872b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 873b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 8742ca6e920SHong Zhang if (ts->vecs_sensip) { 875b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 8762ca6e920SHong Zhang } 877b5b4017aSHong Zhang if (ts->vecs_sensi2) { 878b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 879b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 880b5b4017aSHong Zhang } 881b5b4017aSHong Zhang if (ts->vecs_sensip2) { 882b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 883b5b4017aSHong Zhang } 884316643e7SJed Brown PetscFunctionReturn(0); 885316643e7SJed Brown } 886316643e7SJed Brown 8874416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 888316643e7SJed Brown { 889316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 890316643e7SJed Brown PetscErrorCode ierr; 891316643e7SJed Brown 892316643e7SJed Brown PetscFunctionBegin; 893e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 894316643e7SJed Brown { 8950298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 8960298fd71SBarry 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); 89703842d09SLisandro 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); 898316643e7SJed Brown } 899316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 900316643e7SJed Brown PetscFunctionReturn(0); 901316643e7SJed Brown } 902316643e7SJed Brown 903316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 904316643e7SJed Brown { 905316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 906ace3abfcSBarry Smith PetscBool iascii; 907316643e7SJed Brown PetscErrorCode ierr; 908316643e7SJed Brown 909316643e7SJed Brown PetscFunctionBegin; 910251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 911316643e7SJed Brown if (iascii) { 9127c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 913ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 914316643e7SJed Brown } 915316643e7SJed Brown PetscFunctionReturn(0); 916316643e7SJed Brown } 917316643e7SJed Brown 918560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 9190de4c49aSJed Brown { 9200de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 9210de4c49aSJed Brown 9220de4c49aSJed Brown PetscFunctionBegin; 9230de4c49aSJed Brown *theta = th->Theta; 9240de4c49aSJed Brown PetscFunctionReturn(0); 9250de4c49aSJed Brown } 9260de4c49aSJed Brown 927560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 9280de4c49aSJed Brown { 9290de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 9300de4c49aSJed Brown 9310de4c49aSJed Brown PetscFunctionBegin; 9327c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 9330de4c49aSJed Brown th->Theta = theta; 9341566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 9350de4c49aSJed Brown PetscFunctionReturn(0); 9360de4c49aSJed Brown } 937eb284becSJed Brown 938560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 93926f2ff8fSLisandro Dalcin { 94026f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 94126f2ff8fSLisandro Dalcin 94226f2ff8fSLisandro Dalcin PetscFunctionBegin; 94326f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 94426f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 94526f2ff8fSLisandro Dalcin } 94626f2ff8fSLisandro Dalcin 947560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 948eb284becSJed Brown { 949eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 950eb284becSJed Brown 951eb284becSJed Brown PetscFunctionBegin; 952eb284becSJed Brown th->endpoint = flg; 953eb284becSJed Brown PetscFunctionReturn(0); 954eb284becSJed Brown } 9550de4c49aSJed Brown 956f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 957f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 958f9c1d6abSBarry Smith { 959f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 960f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 9613fd8ae06SJed Brown const PetscReal one = 1.0; 962f9c1d6abSBarry Smith 963f9c1d6abSBarry Smith PetscFunctionBegin; 9643fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 965f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 966f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 967f9c1d6abSBarry Smith PetscFunctionReturn(0); 968f9c1d6abSBarry Smith } 969f9c1d6abSBarry Smith #endif 970f9c1d6abSBarry Smith 97142682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 97242682096SHong Zhang { 97342682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 97442682096SHong Zhang 97542682096SHong Zhang PetscFunctionBegin; 9761566a47fSLisandro Dalcin if (ns) *ns = 1; 9771566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 97842682096SHong Zhang PetscFunctionReturn(0); 97942682096SHong Zhang } 980f9c1d6abSBarry Smith 981316643e7SJed Brown /* ------------------------------------------------------------ */ 982316643e7SJed Brown /*MC 98396f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 984316643e7SJed Brown 985316643e7SJed Brown Level: beginner 986316643e7SJed Brown 9874eb428fdSBarry Smith Options Database: 988c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 989c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 99003842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 9914eb428fdSBarry Smith 992eb284becSJed Brown Notes: 9930c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 9940c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 9954eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 9964eb428fdSBarry Smith 997eb284becSJed Brown This method can be applied to DAE. 998eb284becSJed Brown 999eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 1000eb284becSJed Brown 1001eb284becSJed Brown .vb 1002eb284becSJed Brown Theta | Theta 1003eb284becSJed Brown ------------- 1004eb284becSJed Brown | 1 1005eb284becSJed Brown .ve 1006eb284becSJed Brown 1007eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1008eb284becSJed Brown 1009eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1010eb284becSJed Brown 1011eb284becSJed Brown .vb 1012eb284becSJed Brown 0 | 0 0 1013eb284becSJed Brown 1 | 1-Theta Theta 1014eb284becSJed Brown ------------------- 1015eb284becSJed Brown | 1-Theta Theta 1016eb284becSJed Brown .ve 1017eb284becSJed Brown 1018eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1019eb284becSJed Brown 1020eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1021eb284becSJed Brown 1022eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1023eb284becSJed Brown 10244eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1025eb284becSJed Brown 1026eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1027316643e7SJed Brown 1028316643e7SJed Brown M*/ 10298cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1030316643e7SJed Brown { 1031316643e7SJed Brown TS_Theta *th; 1032316643e7SJed Brown PetscErrorCode ierr; 1033316643e7SJed Brown 1034316643e7SJed Brown PetscFunctionBegin; 1035277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1036316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1037316643e7SJed Brown ts->ops->view = TSView_Theta; 1038316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 103942f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1040316643e7SJed Brown ts->ops->step = TSStep_Theta; 1041cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 10421566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 104324655328SShri ts->ops->rollback = TSRollBack_Theta; 1044316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 10450f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 10460f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1047f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1048f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1049f9c1d6abSBarry Smith #endif 105042682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 105142f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1052b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1053b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 10542ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 1055715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 1056715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 1057316643e7SJed Brown 1058efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1059efd4aadfSBarry Smith 1060b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1061316643e7SJed Brown ts->data = (void*)th; 1062316643e7SJed Brown 1063715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1064715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1065715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1066b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1067715f1b00SHong Zhang 10686f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1069316643e7SJed Brown th->Theta = 0.5; 10701566a47fSLisandro Dalcin th->order = 2; 1071bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1072bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1073bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1074bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1075316643e7SJed Brown PetscFunctionReturn(0); 1076316643e7SJed Brown } 10770de4c49aSJed Brown 10780de4c49aSJed Brown /*@ 10790de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 10800de4c49aSJed Brown 10810de4c49aSJed Brown Not Collective 10820de4c49aSJed Brown 10830de4c49aSJed Brown Input Parameter: 10840de4c49aSJed Brown . ts - timestepping context 10850de4c49aSJed Brown 10860de4c49aSJed Brown Output Parameter: 10870de4c49aSJed Brown . theta - stage abscissa 10880de4c49aSJed Brown 10890de4c49aSJed Brown Note: 10900de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 10910de4c49aSJed Brown 10920de4c49aSJed Brown Level: Advanced 10930de4c49aSJed Brown 10940de4c49aSJed Brown .seealso: TSThetaSetTheta() 10950de4c49aSJed Brown @*/ 10967087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 10970de4c49aSJed Brown { 10984ac538c5SBarry Smith PetscErrorCode ierr; 10990de4c49aSJed Brown 11000de4c49aSJed Brown PetscFunctionBegin; 1101afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 11020de4c49aSJed Brown PetscValidPointer(theta,2); 11034ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 11040de4c49aSJed Brown PetscFunctionReturn(0); 11050de4c49aSJed Brown } 11060de4c49aSJed Brown 11070de4c49aSJed Brown /*@ 11080de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 11090de4c49aSJed Brown 11100de4c49aSJed Brown Not Collective 11110de4c49aSJed Brown 11120de4c49aSJed Brown Input Parameter: 11130de4c49aSJed Brown + ts - timestepping context 11140de4c49aSJed Brown - theta - stage abscissa 11150de4c49aSJed Brown 11160de4c49aSJed Brown Options Database: 11170de4c49aSJed Brown . -ts_theta_theta <theta> 11180de4c49aSJed Brown 11190de4c49aSJed Brown Level: Intermediate 11200de4c49aSJed Brown 11210de4c49aSJed Brown .seealso: TSThetaGetTheta() 11220de4c49aSJed Brown @*/ 11237087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 11240de4c49aSJed Brown { 11254ac538c5SBarry Smith PetscErrorCode ierr; 11260de4c49aSJed Brown 11270de4c49aSJed Brown PetscFunctionBegin; 1128afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 11294ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 11300de4c49aSJed Brown PetscFunctionReturn(0); 11310de4c49aSJed Brown } 1132f33bbcb6SJed Brown 113326f2ff8fSLisandro Dalcin /*@ 113426f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 113526f2ff8fSLisandro Dalcin 113626f2ff8fSLisandro Dalcin Not Collective 113726f2ff8fSLisandro Dalcin 113826f2ff8fSLisandro Dalcin Input Parameter: 113926f2ff8fSLisandro Dalcin . ts - timestepping context 114026f2ff8fSLisandro Dalcin 114126f2ff8fSLisandro Dalcin Output Parameter: 114226f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 114326f2ff8fSLisandro Dalcin 114426f2ff8fSLisandro Dalcin Level: Advanced 114526f2ff8fSLisandro Dalcin 114626f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 114726f2ff8fSLisandro Dalcin @*/ 114826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 114926f2ff8fSLisandro Dalcin { 115026f2ff8fSLisandro Dalcin PetscErrorCode ierr; 115126f2ff8fSLisandro Dalcin 115226f2ff8fSLisandro Dalcin PetscFunctionBegin; 115326f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 115426f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1155163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 115626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 115726f2ff8fSLisandro Dalcin } 115826f2ff8fSLisandro Dalcin 1159eb284becSJed Brown /*@ 1160eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1161eb284becSJed Brown 1162eb284becSJed Brown Not Collective 1163eb284becSJed Brown 1164eb284becSJed Brown Input Parameter: 1165eb284becSJed Brown + ts - timestepping context 1166eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1167eb284becSJed Brown 1168eb284becSJed Brown Options Database: 1169eb284becSJed Brown . -ts_theta_endpoint <flg> 1170eb284becSJed Brown 1171eb284becSJed Brown Level: Intermediate 1172eb284becSJed Brown 1173eb284becSJed Brown .seealso: TSTHETA, TSCN 1174eb284becSJed Brown @*/ 1175eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1176eb284becSJed Brown { 1177eb284becSJed Brown PetscErrorCode ierr; 1178eb284becSJed Brown 1179eb284becSJed Brown PetscFunctionBegin; 1180eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1181eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1182eb284becSJed Brown PetscFunctionReturn(0); 1183eb284becSJed Brown } 1184eb284becSJed Brown 1185f33bbcb6SJed Brown /* 1186f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1187f33bbcb6SJed Brown * The creation functions for these specializations are below. 1188f33bbcb6SJed Brown */ 1189f33bbcb6SJed Brown 11901566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 11911566a47fSLisandro Dalcin { 11921566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 11931566a47fSLisandro Dalcin PetscErrorCode ierr; 11941566a47fSLisandro Dalcin 11951566a47fSLisandro Dalcin PetscFunctionBegin; 11961566a47fSLisandro 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"); 11971566a47fSLisandro 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"); 11981566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 11991566a47fSLisandro Dalcin PetscFunctionReturn(0); 12001566a47fSLisandro Dalcin } 12011566a47fSLisandro Dalcin 1202f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1203f33bbcb6SJed Brown { 1204f33bbcb6SJed Brown PetscFunctionBegin; 1205f33bbcb6SJed Brown PetscFunctionReturn(0); 1206f33bbcb6SJed Brown } 1207f33bbcb6SJed Brown 1208f33bbcb6SJed Brown /*MC 1209f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1210f33bbcb6SJed Brown 1211f33bbcb6SJed Brown Level: beginner 1212f33bbcb6SJed Brown 12134eb428fdSBarry Smith Notes: 1214c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 12154eb428fdSBarry Smith 12161566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 12174eb428fdSBarry Smith 1218f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1219f33bbcb6SJed Brown 1220f33bbcb6SJed Brown M*/ 12218cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1222f33bbcb6SJed Brown { 1223f33bbcb6SJed Brown PetscErrorCode ierr; 1224f33bbcb6SJed Brown 1225f33bbcb6SJed Brown PetscFunctionBegin; 1226f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1227f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 12281566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 12290c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1230f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1231f33bbcb6SJed Brown PetscFunctionReturn(0); 1232f33bbcb6SJed Brown } 1233f33bbcb6SJed Brown 12341566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 12351566a47fSLisandro Dalcin { 12361566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 12371566a47fSLisandro Dalcin PetscErrorCode ierr; 12381566a47fSLisandro Dalcin 12391566a47fSLisandro Dalcin PetscFunctionBegin; 12401566a47fSLisandro 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"); 12411566a47fSLisandro 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"); 12421566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 12431566a47fSLisandro Dalcin PetscFunctionReturn(0); 12441566a47fSLisandro Dalcin } 12451566a47fSLisandro Dalcin 1246f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1247f33bbcb6SJed Brown { 1248f33bbcb6SJed Brown PetscFunctionBegin; 1249f33bbcb6SJed Brown PetscFunctionReturn(0); 1250f33bbcb6SJed Brown } 1251f33bbcb6SJed Brown 1252f33bbcb6SJed Brown /*MC 1253f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1254f33bbcb6SJed Brown 1255f33bbcb6SJed Brown Level: beginner 1256f33bbcb6SJed Brown 1257f33bbcb6SJed Brown Notes: 12587cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 12597cf5af47SJed Brown 12607cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1261f33bbcb6SJed Brown 1262f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1263f33bbcb6SJed Brown 1264f33bbcb6SJed Brown M*/ 12658cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1266f33bbcb6SJed Brown { 1267f33bbcb6SJed Brown PetscErrorCode ierr; 1268f33bbcb6SJed Brown 1269f33bbcb6SJed Brown PetscFunctionBegin; 1270f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1271f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1272eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 12730c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1274f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1275f33bbcb6SJed Brown PetscFunctionReturn(0); 1276f33bbcb6SJed Brown } 1277