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; 133*cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 134022c081aSHong Zhang PetscErrorCode ierr; 135022c081aSHong Zhang 136022c081aSHong Zhang PetscFunctionBegin; 137022c081aSHong Zhang if (th->endpoint) { 138022c081aSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 139022c081aSHong Zhang if (th->Theta!=1.0) { 140*cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 141*cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 142022c081aSHong Zhang } 143*cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 144*cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 145022c081aSHong Zhang } else { 146*cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 147*cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 148022c081aSHong Zhang } 149022c081aSHong Zhang PetscFunctionReturn(0); 150022c081aSHong Zhang } 151022c081aSHong Zhang 152b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 153b1cb36f3SHong Zhang { 154b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 155*cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 156b1cb36f3SHong Zhang PetscErrorCode ierr; 157b1cb36f3SHong Zhang 158b1cb36f3SHong Zhang PetscFunctionBegin; 159b1cb36f3SHong Zhang /* backup cost integral */ 160*cd4cee2dSHong Zhang ierr = VecCopy(quadts->vec_sol,th->VecCostIntegral0);CHKERRQ(ierr); 161022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 162b1cb36f3SHong Zhang PetscFunctionReturn(0); 163b1cb36f3SHong Zhang } 164b1cb36f3SHong Zhang 165b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 166b1cb36f3SHong Zhang { 167b1cb36f3SHong Zhang PetscErrorCode ierr; 168b1cb36f3SHong Zhang 169b1cb36f3SHong Zhang PetscFunctionBegin; 170022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 17124655328SShri PetscFunctionReturn(0); 17224655328SShri } 17324655328SShri 174470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x) 1751566a47fSLisandro Dalcin { 1761566a47fSLisandro Dalcin PetscInt nits,lits; 1771566a47fSLisandro Dalcin PetscErrorCode ierr; 1781566a47fSLisandro Dalcin 1791566a47fSLisandro Dalcin PetscFunctionBegin; 1801566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1811566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1821566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1831566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1841566a47fSLisandro Dalcin PetscFunctionReturn(0); 1851566a47fSLisandro Dalcin } 1861566a47fSLisandro Dalcin 187193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 188316643e7SJed Brown { 189316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1901566a47fSLisandro Dalcin PetscInt rejections = 0; 1914957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 1921566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 193051f2191SLisandro Dalcin PetscErrorCode ierr; 194316643e7SJed Brown 195316643e7SJed Brown PetscFunctionBegin; 1961566a47fSLisandro Dalcin if (!ts->steprollback) { 1972ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 1983b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 1991566a47fSLisandro Dalcin } 2001566a47fSLisandro Dalcin 2011566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 2021566a47fSLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 203715f1b00SHong Zhang 2041566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 2051566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 206316643e7SJed Brown 207be5899b3SLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 208fecfb714SLisandro Dalcin if (th->extrapolate && !ts->steprestart) { 209be5899b3SLisandro Dalcin ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 210be5899b3SLisandro Dalcin } 211eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 212eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2131566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2141566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2151566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2161566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2171566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 218eb284becSJed Brown } 219be5899b3SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 220470880abSPatrick Sanan ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 221be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 2221566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 223fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 224051f2191SLisandro Dalcin 225051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2261566a47fSLisandro Dalcin if (th->endpoint) { 2271566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2281566a47fSLisandro Dalcin } else { 229cb3a5882SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 2301566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2311566a47fSLisandro Dalcin } 2321566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2331566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2341566a47fSLisandro Dalcin if (!accept) { 2351566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 236be5899b3SLisandro Dalcin ts->time_step = next_time_step; 237be5899b3SLisandro Dalcin goto reject_step; 238051f2191SLisandro Dalcin } 2393b1890cdSShri Abhyankar 240715f1b00SHong Zhang if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 241b1cb36f3SHong Zhang th->ptime = ts->ptime; 242b1cb36f3SHong Zhang th->time_step = ts->time_step; 24317145e75SHong Zhang } 2441566a47fSLisandro Dalcin 2452b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 246cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 247051f2191SLisandro Dalcin break; 248051f2191SLisandro Dalcin 249051f2191SLisandro Dalcin reject_step: 250fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 2511566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 252051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2531566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2543b1890cdSShri Abhyankar } 2553b1890cdSShri Abhyankar } 256316643e7SJed Brown PetscFunctionReturn(0); 257316643e7SJed Brown } 258316643e7SJed Brown 259c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 260c9681e5dSHong Zhang { 261c9681e5dSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 262*cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 263c9681e5dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 264c9681e5dSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 265c9681e5dSHong Zhang PetscInt nadj; 266*cd4cee2dSHong Zhang Mat J,Jpre,quadJ,quadJp; 267c9681e5dSHong Zhang KSP ksp; 268c9681e5dSHong Zhang PetscReal shift; 269c9681e5dSHong Zhang PetscScalar *xarr; 270c9681e5dSHong Zhang PetscErrorCode ierr; 271c9681e5dSHong Zhang 272c9681e5dSHong Zhang PetscFunctionBegin; 273c9681e5dSHong Zhang th->status = TS_STEP_INCOMPLETE; 274c9681e5dSHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 275*cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 276*cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 277*cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 278c9681e5dSHong Zhang 279c9681e5dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 280c9681e5dSHong Zhang th->stage_time = ts->ptime; /* time_step is negative*/ 281c9681e5dSHong Zhang th->ptime = ts->ptime + ts->time_step; 282c9681e5dSHong Zhang th->time_step = -ts->time_step; 283c9681e5dSHong Zhang 28433f72c5dSHong Zhang /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */ 285*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 286*cd4cee2dSHong Zhang 287c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 288c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 289c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */ 290*cd4cee2dSHong Zhang if (quadJ) { 291*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 292*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 293*cd4cee2dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr); 294*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 295*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 296c9681e5dSHong Zhang } 297c9681e5dSHong Zhang } 298c9681e5dSHong Zhang 299c9681e5dSHong Zhang /* Build LHS for first-order adjoint */ 300c9681e5dSHong Zhang shift = 1./th->time_step; 301*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 302*cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 303c9681e5dSHong Zhang 304c9681e5dSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 305c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 306c9681e5dSHong Zhang KSPConvergedReason kspreason; 307c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 308c9681e5dSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 309c9681e5dSHong Zhang if (kspreason < 0) { 310c9681e5dSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 311c9681e5dSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 312c9681e5dSHong Zhang } 313c9681e5dSHong Zhang } 314c9681e5dSHong Zhang 315c9681e5dSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 316c9681e5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 317c9681e5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 318c9681e5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 319c9681e5dSHong Zhang /* lambda_s^T F_UU w_1 */ 320c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 321c9681e5dSHong Zhang /* lambda_s^T F_UP w_2 */ 322c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 323c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 324c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 325c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 32633f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 327c9681e5dSHong Zhang if (ts->vecs_fup) { 32833f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 329c9681e5dSHong Zhang } 330c9681e5dSHong Zhang } 331c9681e5dSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 332c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 33305c9950eSHong Zhang KSPConvergedReason kspreason; 334c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 33505c9950eSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 33605c9950eSHong Zhang if (kspreason < 0) { 33705c9950eSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 33805c9950eSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 33905c9950eSHong Zhang } 340c9681e5dSHong Zhang } 341c9681e5dSHong Zhang } 342c9681e5dSHong Zhang 343c9681e5dSHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 344c9681e5dSHong Zhang shift = 0.0; 345*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_U */ 346c9681e5dSHong Zhang ierr = MatScale(J,-1.);CHKERRQ(ierr); 347c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 34833f72c5dSHong Zhang /* Add f_U \lambda_s to the original RHS */ 349c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 350c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr); 351c9681e5dSHong Zhang ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 352c9681e5dSHong Zhang if (ts->vecs_sensi2) { 353c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 354c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr); 355c9681e5dSHong Zhang ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 356c9681e5dSHong Zhang } 357c9681e5dSHong Zhang } 358c9681e5dSHong Zhang if (ts->vecs_sensip) { 359c9681e5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */ 360*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 361c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 362c9681e5dSHong Zhang /* lambda_s^T F_PU w_1 */ 363c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 36433f72c5dSHong Zhang /* lambda_s^T F_PP w_2 */ 365c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 366c9681e5dSHong Zhang } 367*cd4cee2dSHong Zhang 368c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 369c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 370c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 371*cd4cee2dSHong Zhang if (quadJp) { 372*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 373*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 374*cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr); 375*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 376*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 37733f72c5dSHong Zhang } 378c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 379c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 380c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 381c9681e5dSHong Zhang if (ts->vecs_fpu) { 382c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 383c9681e5dSHong Zhang } 384c9681e5dSHong Zhang if (ts->vecs_fpp) { 385c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 386c9681e5dSHong Zhang } 387c9681e5dSHong Zhang } 388c9681e5dSHong Zhang } 389c9681e5dSHong Zhang } 390c9681e5dSHong Zhang 391c9681e5dSHong Zhang if (ts->vecs_sensi2) { 392c9681e5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 393c9681e5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 394c9681e5dSHong Zhang } 395c9681e5dSHong Zhang th->status = TS_STEP_COMPLETE; 396c9681e5dSHong Zhang PetscFunctionReturn(0); 397c9681e5dSHong Zhang } 398c9681e5dSHong Zhang 39942f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 4002ca6e920SHong Zhang { 4012ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 402*cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 403b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 404b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 4052ca6e920SHong Zhang PetscInt nadj; 406*cd4cee2dSHong Zhang Mat J,Jpre,quadJ,quadJp; 4072ca6e920SHong Zhang KSP ksp; 4082ca6e920SHong Zhang PetscReal shift; 409b5b4017aSHong Zhang PetscScalar *xarr; 410b5b4017aSHong Zhang PetscErrorCode ierr; 4112ca6e920SHong Zhang 4122ca6e920SHong Zhang PetscFunctionBegin; 413c9681e5dSHong Zhang if (th->Theta == 1.) { 414c9681e5dSHong Zhang ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr); 415c9681e5dSHong Zhang PetscFunctionReturn(0); 416c9681e5dSHong Zhang } 4172ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 418302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 419*cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 420*cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 421*cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 422b5e0532dSHong Zhang 423b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 424022c081aSHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 425b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 426022c081aSHong Zhang th->time_step = -ts->time_step; 4272ca6e920SHong Zhang 428b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 42933f72c5dSHong Zhang /* Cost function has an integral term */ 43005755b9cSHong Zhang if (th->endpoint) { 431*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr); 43205755b9cSHong Zhang } else { 433*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 43405755b9cSHong Zhang } 43533f72c5dSHong Zhang 436abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 437b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 438022c081aSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr); 439*cd4cee2dSHong Zhang if (quadJ) { 440*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 441*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 442*cd4cee2dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr); 443*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 444*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 44536eaed60SHong Zhang } 4462ca6e920SHong Zhang } 4473c54f92cSHong Zhang 448b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 449022c081aSHong Zhang shift = 1./(th->Theta*th->time_step); 4503c54f92cSHong Zhang if (th->endpoint) { 451*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 4523c54f92cSHong Zhang } else { 453*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 4543c54f92cSHong Zhang } 455*cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 4562ca6e920SHong Zhang 457b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 458abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 4591d14d03bSHong Zhang KSPConvergedReason kspreason; 460b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 4611d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 4621d14d03bSHong Zhang if (kspreason < 0) { 4631d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 4641d14d03bSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 4651d14d03bSHong Zhang } 4662ca6e920SHong Zhang } 4673c54f92cSHong Zhang 4681d14d03bSHong Zhang /* Second-order adjoint */ 469b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 470b5b4017aSHong Zhang if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 471b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 472b5b4017aSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 473b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 474b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 475b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 476b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 477b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 478b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 479b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 480b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 481b5b4017aSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 4821d14d03bSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 48333f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 484b5b4017aSHong Zhang if (ts->vecs_fup) { 48533f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 486b5b4017aSHong Zhang } 487b5b4017aSHong Zhang } 488b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 489b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 4901d14d03bSHong Zhang KSPConvergedReason kspreason; 4911d14d03bSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 4921d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 4931d14d03bSHong Zhang if (kspreason < 0) { 4941d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 4951d14d03bSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 4961d14d03bSHong Zhang } 497b5b4017aSHong Zhang } 498b5b4017aSHong Zhang } 499b5b4017aSHong Zhang 50036eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 5011d14d03bSHong Zhang if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 502022c081aSHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 503*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 50433f72c5dSHong Zhang /* R_U at t_n */ 505*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr); 506abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 507b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 508*cd4cee2dSHong Zhang if (quadJ) { 509*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 510*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 511*cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr); 512*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 513*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 514b5b4017aSHong Zhang } 5151d14d03bSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 516b5b4017aSHong Zhang } 5171d14d03bSHong Zhang 5181d14d03bSHong Zhang /* Second-order adjoint */ 5191d14d03bSHong Zhang if (ts->vecs_sensi2) { /* U_n */ 520b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 521b5b4017aSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 522b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 523b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 524b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 525b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 52633f72c5dSHong Zhang ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr); 527b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 528b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 529b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 53033f72c5dSHong Zhang /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2 */ 531b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 53233f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 533b5b4017aSHong Zhang if (ts->vecs_fup) { 53433f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 535b5b4017aSHong Zhang } 5361d14d03bSHong Zhang ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr); 537b5b4017aSHong Zhang } 538b5e0532dSHong Zhang } 5393fd52205SHong Zhang 5403fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 541b5b4017aSHong Zhang /* U_{n+1} */ 54233f72c5dSHong Zhang shift = -1./(th->Theta*th->time_step); 54333f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 544*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 545abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 546b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 54733f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 5483fd52205SHong Zhang } 54933f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 55033f72c5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 55133f72c5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 55233f72c5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 553b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 554b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 55533f72c5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 55633f72c5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 55733f72c5dSHong Zhang 558b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 559b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 560b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 56133f72c5dSHong Zhang /* Mu2 <- Mu2 + h theta F_P^T Lambda_s + h theta (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */ 562b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 56333f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 564b5b4017aSHong Zhang if (ts->vecs_fpu) { 56533f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 566b5b4017aSHong Zhang } 567b5b4017aSHong Zhang if (ts->vecs_fpp) { 56833f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 569b5b4017aSHong Zhang } 570b5b4017aSHong Zhang } 571b5b4017aSHong Zhang } 572b5b4017aSHong Zhang 573b5b4017aSHong Zhang /* U_s */ 57433f72c5dSHong Zhang shift = 1./((th->Theta-1.0)*th->time_step); 57533f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 576*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr); 577abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 578b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 57933f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 58033f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 58133f72c5dSHong Zhang /* Get w1 at t_n from TLM matrix */ 58233f72c5dSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 58333f72c5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 584b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 58533f72c5dSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 58633f72c5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 58733f72c5dSHong Zhang ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr); 588b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 58933f72c5dSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 590b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 59133f72c5dSHong Zhang /* Mu2 <- Mu2 + h(1-theta) F_P^T Lambda_s + h(1-theta) (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */ 592b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 59333f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 594b5b4017aSHong Zhang if (ts->vecs_fpu) { 59533f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 596b5e0532dSHong Zhang } 597b5b4017aSHong Zhang if (ts->vecs_fpp) { 59833f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 599b5b4017aSHong Zhang } 60036eaed60SHong Zhang } 60136eaed60SHong Zhang } 6023fd52205SHong Zhang } 603b5e0532dSHong Zhang } 6043fd52205SHong Zhang } else { /* one-stage case */ 6053c54f92cSHong Zhang shift = 0.0; 606*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 607*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 608abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 609b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 610022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 611*cd4cee2dSHong Zhang if (quadJ) { 612*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 613*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 614*cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vec_drdu_col);CHKERRQ(ierr); 615*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 616*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 61736eaed60SHong Zhang } 6182ca6e920SHong Zhang } 6193fd52205SHong Zhang if (ts->vecs_sensip) { 62033f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 621*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 622abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 62333f72c5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 62433f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 625*cd4cee2dSHong Zhang if (quadJp) { 626*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 627*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 628*cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr); 629*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 630*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 63136eaed60SHong Zhang } 63236eaed60SHong Zhang } 6333fd52205SHong Zhang } 6342ca6e920SHong Zhang } 6352ca6e920SHong Zhang 6362ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6372ca6e920SHong Zhang PetscFunctionReturn(0); 6382ca6e920SHong Zhang } 6392ca6e920SHong Zhang 640cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 641cd652676SJed Brown { 642cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 643be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 644cd652676SJed Brown PetscErrorCode ierr; 645cd652676SJed Brown 646cd652676SJed Brown PetscFunctionBegin; 647a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 648be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 649be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 650cd652676SJed Brown PetscFunctionReturn(0); 651cd652676SJed Brown } 652cd652676SJed Brown 6531566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 6541566a47fSLisandro Dalcin { 6551566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 6561566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 6571566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 6587453f775SEmil Constantinescu PetscReal wltea,wlter; 6591566a47fSLisandro Dalcin PetscErrorCode ierr; 6601566a47fSLisandro Dalcin 6611566a47fSLisandro Dalcin PetscFunctionBegin; 6622ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 6631566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 664fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 6651566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 666fecfb714SLisandro Dalcin { 667be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 668be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 6691566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 6701566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 6711566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 6721566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 6731566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 6747453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 6751566a47fSLisandro Dalcin } 6761566a47fSLisandro Dalcin if (order) *order = 2; 6771566a47fSLisandro Dalcin PetscFunctionReturn(0); 6781566a47fSLisandro Dalcin } 6791566a47fSLisandro Dalcin 6801566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 6811566a47fSLisandro Dalcin { 6821566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 68313b1b0ffSHong Zhang PetscInt ncost; 6841566a47fSLisandro Dalcin PetscErrorCode ierr; 6851566a47fSLisandro Dalcin 6861566a47fSLisandro Dalcin PetscFunctionBegin; 6871566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 688*cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 689*cd4cee2dSHong Zhang ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr); 6901566a47fSLisandro Dalcin } 691715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 69213b1b0ffSHong Zhang if (ts->mat_sensip) { 69313b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 6946f92202bSHong Zhang } 6956f92202bSHong Zhang if (ts->vecs_integral_sensip) { 6966f92202bSHong Zhang for (ncost=0;ncost<ts->numcost;ncost++) { 6976f92202bSHong Zhang ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 6986f92202bSHong Zhang } 6996f92202bSHong Zhang } 700715f1b00SHong Zhang PetscFunctionReturn(0); 701715f1b00SHong Zhang } 702715f1b00SHong Zhang 703715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 704715f1b00SHong Zhang { 705715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 706*cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 70713b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 708b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 70913b1b0ffSHong Zhang PetscInt ncost,ntlm; 710715f1b00SHong Zhang KSP ksp; 711*cd4cee2dSHong Zhang Mat J,Jpre,quadJ,quadJp; 712715f1b00SHong Zhang PetscReal shift; 71313b1b0ffSHong Zhang PetscScalar *barr,*xarr; 714715f1b00SHong Zhang PetscErrorCode ierr; 715715f1b00SHong Zhang 716715f1b00SHong Zhang PetscFunctionBegin; 71713b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 71813b1b0ffSHong Zhang 7196f92202bSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 7206f92202bSHong Zhang if (ts->vecs_integral_sensip) { 7216f92202bSHong Zhang ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 7226f92202bSHong Zhang } 7236f92202bSHong Zhang } 7246f92202bSHong Zhang 725715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 726*cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 727*cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 728*cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 729715f1b00SHong Zhang 730715f1b00SHong Zhang /* Build RHS */ 731715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 732715f1b00SHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 733*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 73413b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 73513b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 736715f1b00SHong Zhang 737022c081aSHong Zhang /* Add the f_p forcing terms */ 7380e11c21fSHong Zhang if (ts->Jacp) { 73933f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 74033f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 74133f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 74233f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 7430e11c21fSHong Zhang } 744715f1b00SHong Zhang } else { /* 1-stage method */ 745715f1b00SHong Zhang shift = 0.0; 746*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 74713b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 74813b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 74913b1b0ffSHong Zhang 750022c081aSHong Zhang /* Add the f_p forcing terms */ 7510e11c21fSHong Zhang if (ts->Jacp) { 75233f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 75333f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 754715f1b00SHong Zhang } 7550e11c21fSHong Zhang } 756715f1b00SHong Zhang 757715f1b00SHong Zhang /* Build LHS */ 758715f1b00SHong Zhang shift = 1/(th->Theta*th->time_step); 759715f1b00SHong Zhang if (th->endpoint) { 760*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 761715f1b00SHong Zhang } else { 762*cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 763715f1b00SHong Zhang } 764*cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 765715f1b00SHong Zhang 766715f1b00SHong Zhang /* 767715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 768c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 769715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 770715f1b00SHong Zhang */ 771715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 772715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 773*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr); 774*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr); 775715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 776*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,ncost,&xarr);CHKERRQ(ierr); 777*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 778*cd4cee2dSHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vec_drdu_col,th->VecIntegralSensipTemp);CHKERRQ(ierr); 779*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 780*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 781*cd4cee2dSHong Zhang if (quadJp) { 782*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,ncost,&xarr);CHKERRQ(ierr); 783*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 784*cd4cee2dSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vec_drdp_col);CHKERRQ(ierr); 785*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 786*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 7870e11c21fSHong Zhang } 788715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 789715f1b00SHong Zhang } 790715f1b00SHong Zhang } 791715f1b00SHong Zhang } 792715f1b00SHong Zhang 793715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 794022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 795b5b4017aSHong Zhang KSPConvergedReason kspreason; 79613b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 797b5b4017aSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 798715f1b00SHong Zhang if (th->endpoint) { 79913b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 800b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 801b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 802b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 80313b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 804715f1b00SHong Zhang } else { 805b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 806715f1b00SHong Zhang } 807b5b4017aSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 808b5b4017aSHong Zhang if (kspreason < 0) { 809b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 810b5b4017aSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 811b5b4017aSHong Zhang } 812b5b4017aSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 81313b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 814715f1b00SHong Zhang } 815715f1b00SHong Zhang 816b5b4017aSHong Zhang 81713b1b0ffSHong Zhang /* 81813b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 819c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 82013b1b0ffSHong Zhang */ 82113b1b0ffSHong Zhang if (ts->vecs_integral_sensip) { 82213b1b0ffSHong Zhang if (!th->endpoint) { 82313b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 824*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 825*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 82613b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 827*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,ncost,&xarr);CHKERRQ(ierr); 828*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 829*cd4cee2dSHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vec_drdu_col,th->VecIntegralSensipTemp);CHKERRQ(ierr); 830*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 831*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 832*cd4cee2dSHong Zhang if (quadJp) { 833*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,ncost,&xarr);CHKERRQ(ierr); 834*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 835*cd4cee2dSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vec_drdp_col);CHKERRQ(ierr); 836*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 837*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 8380e11c21fSHong Zhang } 83913b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 84013b1b0ffSHong Zhang } 84113b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 84213b1b0ffSHong Zhang } else { 843*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr); 844*cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 84513b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 846*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,ncost,&xarr);CHKERRQ(ierr); 847*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 848*cd4cee2dSHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vec_drdu_col,th->VecIntegralSensipTemp);CHKERRQ(ierr); 849*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 850*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 851*cd4cee2dSHong Zhang if (quadJp) { 852*cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,ncost,&xarr);CHKERRQ(ierr); 853*cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 854*cd4cee2dSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vec_drdp_col);CHKERRQ(ierr); 855*cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 856*cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 8570e11c21fSHong Zhang } 85813b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 85913b1b0ffSHong Zhang } 86013b1b0ffSHong Zhang } 86113b1b0ffSHong Zhang } else { 86213b1b0ffSHong Zhang if (!th->endpoint) { 86313b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 864715f1b00SHong Zhang } 865715f1b00SHong Zhang } 8661566a47fSLisandro Dalcin PetscFunctionReturn(0); 8671566a47fSLisandro Dalcin } 8681566a47fSLisandro Dalcin 8696affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 8706affb6f8SHong Zhang { 8716affb6f8SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 8726affb6f8SHong Zhang 8736affb6f8SHong Zhang PetscFunctionBegin; 8746affb6f8SHong Zhang if (ns) *ns = 1; 8756affb6f8SHong Zhang if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 8766affb6f8SHong Zhang PetscFunctionReturn(0); 8776affb6f8SHong Zhang } 8786affb6f8SHong Zhang 879316643e7SJed Brown /*------------------------------------------------------------*/ 880277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 881316643e7SJed Brown { 882316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 883316643e7SJed Brown PetscErrorCode ierr; 884316643e7SJed Brown 885316643e7SJed Brown PetscFunctionBegin; 8866bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 8876bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 8883b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 889eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 8901566a47fSLisandro Dalcin 8911566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 8921566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 8931566a47fSLisandro Dalcin 894b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 895277b19d0SLisandro Dalcin PetscFunctionReturn(0); 896277b19d0SLisandro Dalcin } 897277b19d0SLisandro Dalcin 898ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts) 899ece46509SHong Zhang { 900ece46509SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 901ece46509SHong Zhang PetscErrorCode ierr; 902ece46509SHong Zhang 903ece46509SHong Zhang PetscFunctionBegin; 904ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 905ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 906ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 907ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 908ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 909ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 910ece46509SHong Zhang PetscFunctionReturn(0); 911ece46509SHong Zhang } 912ece46509SHong Zhang 913277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 914277b19d0SLisandro Dalcin { 915277b19d0SLisandro Dalcin PetscErrorCode ierr; 916277b19d0SLisandro Dalcin 917277b19d0SLisandro Dalcin PetscFunctionBegin; 918277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 919b3a6b972SJed Brown if (ts->dm) { 920b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 921b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 922b3a6b972SJed Brown } 923277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 924bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 925bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 926bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 927bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 928316643e7SJed Brown PetscFunctionReturn(0); 929316643e7SJed Brown } 930316643e7SJed Brown 931316643e7SJed Brown /* 932316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9332b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 934316643e7SJed Brown */ 9350f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 936316643e7SJed Brown { 937316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 938316643e7SJed Brown PetscErrorCode ierr; 9397445fe48SJed Brown Vec X0,Xdot; 9407445fe48SJed Brown DM dm,dmsave; 9411566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 942316643e7SJed Brown 943316643e7SJed Brown PetscFunctionBegin; 9447445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 9455a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9467445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 947b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 9487445fe48SJed Brown 9497445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9507445fe48SJed Brown dmsave = ts->dm; 9517445fe48SJed Brown ts->dm = dm; 9527445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 9537445fe48SJed Brown ts->dm = dmsave; 9540d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 955316643e7SJed Brown PetscFunctionReturn(0); 956316643e7SJed Brown } 957316643e7SJed Brown 958d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 959316643e7SJed Brown { 960316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 961316643e7SJed Brown PetscErrorCode ierr; 9627445fe48SJed Brown Vec Xdot; 9637445fe48SJed Brown DM dm,dmsave; 9641566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 965316643e7SJed Brown 966316643e7SJed Brown PetscFunctionBegin; 9677445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 968be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 9690298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 9707445fe48SJed Brown 9717445fe48SJed Brown dmsave = ts->dm; 9727445fe48SJed Brown ts->dm = dm; 973d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 9747445fe48SJed Brown ts->dm = dmsave; 9750298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 976316643e7SJed Brown PetscFunctionReturn(0); 977316643e7SJed Brown } 978316643e7SJed Brown 979715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 980715f1b00SHong Zhang { 981715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 982715f1b00SHong Zhang PetscErrorCode ierr; 983715f1b00SHong Zhang 984715f1b00SHong Zhang PetscFunctionBegin; 985022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 98613b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 98713b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 988715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 989715f1b00SHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 990715f1b00SHong Zhang } 9916f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 99213b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 99313b1b0ffSHong Zhang 9946f92202bSHong Zhang if (ts->vecs_integral_sensip) { 9956f92202bSHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 9966f92202bSHong Zhang } 997b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 998715f1b00SHong Zhang PetscFunctionReturn(0); 999715f1b00SHong Zhang } 1000715f1b00SHong Zhang 10017adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts) 10027adebddeSHong Zhang { 10037adebddeSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 10047adebddeSHong Zhang PetscErrorCode ierr; 10057adebddeSHong Zhang 10067adebddeSHong Zhang PetscFunctionBegin; 10077adebddeSHong Zhang if (ts->vecs_integral_sensip) { 10087adebddeSHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 10097adebddeSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 10107adebddeSHong Zhang } 10117adebddeSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 10127adebddeSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 10137adebddeSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 10147adebddeSHong Zhang PetscFunctionReturn(0); 10157adebddeSHong Zhang } 10167adebddeSHong Zhang 1017316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 1018316643e7SJed Brown { 1019316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1020*cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 10212ffb9264SLisandro Dalcin PetscBool match; 1022316643e7SJed Brown PetscErrorCode ierr; 1023316643e7SJed Brown 1024316643e7SJed Brown PetscFunctionBegin; 1025*cd4cee2dSHong Zhang if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */ 1026*cd4cee2dSHong Zhang ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr); 1027b1cb36f3SHong Zhang } 102839d13666SHong Zhang if (!th->X) { 1029316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 103039d13666SHong Zhang } 103139d13666SHong Zhang if (!th->Xdot) { 1032316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 103339d13666SHong Zhang } 103439d13666SHong Zhang if (!th->X0) { 10353b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 103639d13666SHong Zhang } 10371566a47fSLisandro Dalcin if (th->endpoint) { 10381566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 10397445fe48SJed Brown } 10403b1890cdSShri Abhyankar 10411566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 10421566a47fSLisandro Dalcin 10431566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 10441566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10451566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10461566a47fSLisandro Dalcin 10471566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 10481566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 10492ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 10502ffb9264SLisandro Dalcin if (!match) { 10511566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 10521566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 10533b1890cdSShri Abhyankar } 10541566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 1055b4dd3bf9SBarry Smith PetscFunctionReturn(0); 1056b4dd3bf9SBarry Smith } 10570c7376e5SHong Zhang 1058b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1059b4dd3bf9SBarry Smith 1060b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1061b4dd3bf9SBarry Smith { 1062b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 1063b4dd3bf9SBarry Smith PetscErrorCode ierr; 1064b4dd3bf9SBarry Smith 1065b4dd3bf9SBarry Smith PetscFunctionBegin; 1066b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 1067b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 10682ca6e920SHong Zhang if (ts->vecs_sensip) { 1069b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 10702ca6e920SHong Zhang } 1071b5b4017aSHong Zhang if (ts->vecs_sensi2) { 1072b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 1073b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 107467633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 107567633408SHong Zhang if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu; 107667633408SHong Zhang if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup; 1077b5b4017aSHong Zhang } 1078c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 1079c9681e5dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 108067633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 108167633408SHong Zhang if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu; 108267633408SHong Zhang if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp; 1083b5b4017aSHong Zhang } 1084316643e7SJed Brown PetscFunctionReturn(0); 1085316643e7SJed Brown } 1086316643e7SJed Brown 10874416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 1088316643e7SJed Brown { 1089316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1090316643e7SJed Brown PetscErrorCode ierr; 1091316643e7SJed Brown 1092316643e7SJed Brown PetscFunctionBegin; 1093e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 1094316643e7SJed Brown { 10950298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 10960298fd71SBarry 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); 109703842d09SLisandro 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); 1098316643e7SJed Brown } 1099316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1100316643e7SJed Brown PetscFunctionReturn(0); 1101316643e7SJed Brown } 1102316643e7SJed Brown 1103316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 1104316643e7SJed Brown { 1105316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1106ace3abfcSBarry Smith PetscBool iascii; 1107316643e7SJed Brown PetscErrorCode ierr; 1108316643e7SJed Brown 1109316643e7SJed Brown PetscFunctionBegin; 1110251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1111316643e7SJed Brown if (iascii) { 11127c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 1113ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1114316643e7SJed Brown } 1115316643e7SJed Brown PetscFunctionReturn(0); 1116316643e7SJed Brown } 1117316643e7SJed Brown 1118560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 11190de4c49aSJed Brown { 11200de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11210de4c49aSJed Brown 11220de4c49aSJed Brown PetscFunctionBegin; 11230de4c49aSJed Brown *theta = th->Theta; 11240de4c49aSJed Brown PetscFunctionReturn(0); 11250de4c49aSJed Brown } 11260de4c49aSJed Brown 1127560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 11280de4c49aSJed Brown { 11290de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11300de4c49aSJed Brown 11310de4c49aSJed Brown PetscFunctionBegin; 11327c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 11330de4c49aSJed Brown th->Theta = theta; 11341566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 11350de4c49aSJed Brown PetscFunctionReturn(0); 11360de4c49aSJed Brown } 1137eb284becSJed Brown 1138560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 113926f2ff8fSLisandro Dalcin { 114026f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 114126f2ff8fSLisandro Dalcin 114226f2ff8fSLisandro Dalcin PetscFunctionBegin; 114326f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 114426f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 114526f2ff8fSLisandro Dalcin } 114626f2ff8fSLisandro Dalcin 1147560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 1148eb284becSJed Brown { 1149eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1150eb284becSJed Brown 1151eb284becSJed Brown PetscFunctionBegin; 1152eb284becSJed Brown th->endpoint = flg; 1153eb284becSJed Brown PetscFunctionReturn(0); 1154eb284becSJed Brown } 11550de4c49aSJed Brown 1156f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1157f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 1158f9c1d6abSBarry Smith { 1159f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 1160f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 11613fd8ae06SJed Brown const PetscReal one = 1.0; 1162f9c1d6abSBarry Smith 1163f9c1d6abSBarry Smith PetscFunctionBegin; 11643fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 1165f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 1166f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 1167f9c1d6abSBarry Smith PetscFunctionReturn(0); 1168f9c1d6abSBarry Smith } 1169f9c1d6abSBarry Smith #endif 1170f9c1d6abSBarry Smith 117142682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 117242682096SHong Zhang { 117342682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 117442682096SHong Zhang 117542682096SHong Zhang PetscFunctionBegin; 11761566a47fSLisandro Dalcin if (ns) *ns = 1; 11771566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 117842682096SHong Zhang PetscFunctionReturn(0); 117942682096SHong Zhang } 1180f9c1d6abSBarry Smith 1181316643e7SJed Brown /* ------------------------------------------------------------ */ 1182316643e7SJed Brown /*MC 118396f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1184316643e7SJed Brown 1185316643e7SJed Brown Level: beginner 1186316643e7SJed Brown 11874eb428fdSBarry Smith Options Database: 1188c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1189c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 119003842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 11914eb428fdSBarry Smith 1192eb284becSJed Brown Notes: 11930c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 11940c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 11954eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 11964eb428fdSBarry Smith 1197eb284becSJed Brown This method can be applied to DAE. 1198eb284becSJed Brown 1199eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 1200eb284becSJed Brown 1201eb284becSJed Brown .vb 1202eb284becSJed Brown Theta | Theta 1203eb284becSJed Brown ------------- 1204eb284becSJed Brown | 1 1205eb284becSJed Brown .ve 1206eb284becSJed Brown 1207eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1208eb284becSJed Brown 1209eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1210eb284becSJed Brown 1211eb284becSJed Brown .vb 1212eb284becSJed Brown 0 | 0 0 1213eb284becSJed Brown 1 | 1-Theta Theta 1214eb284becSJed Brown ------------------- 1215eb284becSJed Brown | 1-Theta Theta 1216eb284becSJed Brown .ve 1217eb284becSJed Brown 1218eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1219eb284becSJed Brown 1220eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1221eb284becSJed Brown 1222eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1223eb284becSJed Brown 12244eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1225eb284becSJed Brown 1226eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1227316643e7SJed Brown 1228316643e7SJed Brown M*/ 12298cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1230316643e7SJed Brown { 1231316643e7SJed Brown TS_Theta *th; 1232316643e7SJed Brown PetscErrorCode ierr; 1233316643e7SJed Brown 1234316643e7SJed Brown PetscFunctionBegin; 1235277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1236ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1237316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1238316643e7SJed Brown ts->ops->view = TSView_Theta; 1239316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 124042f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1241ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1242316643e7SJed Brown ts->ops->step = TSStep_Theta; 1243cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12441566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 124524655328SShri ts->ops->rollback = TSRollBack_Theta; 1246316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12470f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12480f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1249f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1250f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1251f9c1d6abSBarry Smith #endif 125242682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 125342f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1254b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1255b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12562ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12576affb6f8SHong Zhang 1258715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12597adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1260715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12616affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1262316643e7SJed Brown 1263efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1264efd4aadfSBarry Smith 1265b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1266316643e7SJed Brown ts->data = (void*)th; 1267316643e7SJed Brown 1268715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1269715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1270715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1271b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1272715f1b00SHong Zhang 12736f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1274316643e7SJed Brown th->Theta = 0.5; 12751566a47fSLisandro Dalcin th->order = 2; 1276bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1277bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1278bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1279bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1280316643e7SJed Brown PetscFunctionReturn(0); 1281316643e7SJed Brown } 12820de4c49aSJed Brown 12830de4c49aSJed Brown /*@ 12840de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 12850de4c49aSJed Brown 12860de4c49aSJed Brown Not Collective 12870de4c49aSJed Brown 12880de4c49aSJed Brown Input Parameter: 12890de4c49aSJed Brown . ts - timestepping context 12900de4c49aSJed Brown 12910de4c49aSJed Brown Output Parameter: 12920de4c49aSJed Brown . theta - stage abscissa 12930de4c49aSJed Brown 12940de4c49aSJed Brown Note: 12950de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 12960de4c49aSJed Brown 12970de4c49aSJed Brown Level: Advanced 12980de4c49aSJed Brown 12990de4c49aSJed Brown .seealso: TSThetaSetTheta() 13000de4c49aSJed Brown @*/ 13017087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 13020de4c49aSJed Brown { 13034ac538c5SBarry Smith PetscErrorCode ierr; 13040de4c49aSJed Brown 13050de4c49aSJed Brown PetscFunctionBegin; 1306afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13070de4c49aSJed Brown PetscValidPointer(theta,2); 13084ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 13090de4c49aSJed Brown PetscFunctionReturn(0); 13100de4c49aSJed Brown } 13110de4c49aSJed Brown 13120de4c49aSJed Brown /*@ 13130de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 13140de4c49aSJed Brown 13150de4c49aSJed Brown Not Collective 13160de4c49aSJed Brown 13170de4c49aSJed Brown Input Parameter: 13180de4c49aSJed Brown + ts - timestepping context 13190de4c49aSJed Brown - theta - stage abscissa 13200de4c49aSJed Brown 13210de4c49aSJed Brown Options Database: 13220de4c49aSJed Brown . -ts_theta_theta <theta> 13230de4c49aSJed Brown 13240de4c49aSJed Brown Level: Intermediate 13250de4c49aSJed Brown 13260de4c49aSJed Brown .seealso: TSThetaGetTheta() 13270de4c49aSJed Brown @*/ 13287087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 13290de4c49aSJed Brown { 13304ac538c5SBarry Smith PetscErrorCode ierr; 13310de4c49aSJed Brown 13320de4c49aSJed Brown PetscFunctionBegin; 1333afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13344ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 13350de4c49aSJed Brown PetscFunctionReturn(0); 13360de4c49aSJed Brown } 1337f33bbcb6SJed Brown 133826f2ff8fSLisandro Dalcin /*@ 133926f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 134026f2ff8fSLisandro Dalcin 134126f2ff8fSLisandro Dalcin Not Collective 134226f2ff8fSLisandro Dalcin 134326f2ff8fSLisandro Dalcin Input Parameter: 134426f2ff8fSLisandro Dalcin . ts - timestepping context 134526f2ff8fSLisandro Dalcin 134626f2ff8fSLisandro Dalcin Output Parameter: 134726f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 134826f2ff8fSLisandro Dalcin 134926f2ff8fSLisandro Dalcin Level: Advanced 135026f2ff8fSLisandro Dalcin 135126f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 135226f2ff8fSLisandro Dalcin @*/ 135326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 135426f2ff8fSLisandro Dalcin { 135526f2ff8fSLisandro Dalcin PetscErrorCode ierr; 135626f2ff8fSLisandro Dalcin 135726f2ff8fSLisandro Dalcin PetscFunctionBegin; 135826f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 135926f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1360163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 136126f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 136226f2ff8fSLisandro Dalcin } 136326f2ff8fSLisandro Dalcin 1364eb284becSJed Brown /*@ 1365eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1366eb284becSJed Brown 1367eb284becSJed Brown Not Collective 1368eb284becSJed Brown 1369eb284becSJed Brown Input Parameter: 1370eb284becSJed Brown + ts - timestepping context 1371eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1372eb284becSJed Brown 1373eb284becSJed Brown Options Database: 1374eb284becSJed Brown . -ts_theta_endpoint <flg> 1375eb284becSJed Brown 1376eb284becSJed Brown Level: Intermediate 1377eb284becSJed Brown 1378eb284becSJed Brown .seealso: TSTHETA, TSCN 1379eb284becSJed Brown @*/ 1380eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1381eb284becSJed Brown { 1382eb284becSJed Brown PetscErrorCode ierr; 1383eb284becSJed Brown 1384eb284becSJed Brown PetscFunctionBegin; 1385eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1386eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1387eb284becSJed Brown PetscFunctionReturn(0); 1388eb284becSJed Brown } 1389eb284becSJed Brown 1390f33bbcb6SJed Brown /* 1391f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1392f33bbcb6SJed Brown * The creation functions for these specializations are below. 1393f33bbcb6SJed Brown */ 1394f33bbcb6SJed Brown 13951566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 13961566a47fSLisandro Dalcin { 13971566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 13981566a47fSLisandro Dalcin PetscErrorCode ierr; 13991566a47fSLisandro Dalcin 14001566a47fSLisandro Dalcin PetscFunctionBegin; 14011566a47fSLisandro 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"); 14021566a47fSLisandro 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"); 14031566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14041566a47fSLisandro Dalcin PetscFunctionReturn(0); 14051566a47fSLisandro Dalcin } 14061566a47fSLisandro Dalcin 1407f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1408f33bbcb6SJed Brown { 1409f33bbcb6SJed Brown PetscFunctionBegin; 1410f33bbcb6SJed Brown PetscFunctionReturn(0); 1411f33bbcb6SJed Brown } 1412f33bbcb6SJed Brown 1413f33bbcb6SJed Brown /*MC 1414f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1415f33bbcb6SJed Brown 1416f33bbcb6SJed Brown Level: beginner 1417f33bbcb6SJed Brown 14184eb428fdSBarry Smith Notes: 1419c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 14204eb428fdSBarry Smith 14211566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 14224eb428fdSBarry Smith 1423f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1424f33bbcb6SJed Brown 1425f33bbcb6SJed Brown M*/ 14268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1427f33bbcb6SJed Brown { 1428f33bbcb6SJed Brown PetscErrorCode ierr; 1429f33bbcb6SJed Brown 1430f33bbcb6SJed Brown PetscFunctionBegin; 1431f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1432f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 14331566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 14340c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1435f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1436f33bbcb6SJed Brown PetscFunctionReturn(0); 1437f33bbcb6SJed Brown } 1438f33bbcb6SJed Brown 14391566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 14401566a47fSLisandro Dalcin { 14411566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14421566a47fSLisandro Dalcin PetscErrorCode ierr; 14431566a47fSLisandro Dalcin 14441566a47fSLisandro Dalcin PetscFunctionBegin; 14451566a47fSLisandro 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"); 14461566a47fSLisandro 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"); 14471566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14481566a47fSLisandro Dalcin PetscFunctionReturn(0); 14491566a47fSLisandro Dalcin } 14501566a47fSLisandro Dalcin 1451f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1452f33bbcb6SJed Brown { 1453f33bbcb6SJed Brown PetscFunctionBegin; 1454f33bbcb6SJed Brown PetscFunctionReturn(0); 1455f33bbcb6SJed Brown } 1456f33bbcb6SJed Brown 1457f33bbcb6SJed Brown /*MC 1458f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1459f33bbcb6SJed Brown 1460f33bbcb6SJed Brown Level: beginner 1461f33bbcb6SJed Brown 1462f33bbcb6SJed Brown Notes: 14637cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 14647cf5af47SJed Brown 14657cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1466f33bbcb6SJed Brown 1467f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1468f33bbcb6SJed Brown 1469f33bbcb6SJed Brown M*/ 14708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1471f33bbcb6SJed Brown { 1472f33bbcb6SJed Brown PetscErrorCode ierr; 1473f33bbcb6SJed Brown 1474f33bbcb6SJed Brown PetscFunctionBegin; 1475f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1476f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1477eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 14780c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1479f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1480f33bbcb6SJed Brown PetscFunctionReturn(0); 1481f33bbcb6SJed Brown } 1482