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 257c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 258c9681e5dSHong Zhang { 259c9681e5dSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 260c9681e5dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 261c9681e5dSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 262c9681e5dSHong Zhang PetscInt nadj; 263c9681e5dSHong Zhang Mat J,Jp; 264c9681e5dSHong Zhang KSP ksp; 265c9681e5dSHong Zhang PetscReal shift; 266c9681e5dSHong Zhang PetscScalar *xarr; 267c9681e5dSHong Zhang PetscErrorCode ierr; 268c9681e5dSHong Zhang 269c9681e5dSHong Zhang PetscFunctionBegin; 270c9681e5dSHong Zhang th->status = TS_STEP_INCOMPLETE; 271c9681e5dSHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 272c9681e5dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 273c9681e5dSHong Zhang 274c9681e5dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 275c9681e5dSHong Zhang th->stage_time = ts->ptime; /* time_step is negative*/ 276c9681e5dSHong Zhang th->ptime = ts->ptime + ts->time_step; 277c9681e5dSHong Zhang th->time_step = -ts->time_step; 278c9681e5dSHong Zhang 279c9681e5dSHong Zhang /* Build RHS for first-order adjoint */ 280c9681e5dSHong Zhang if (ts->vec_costintegral) { /* Cost function has an integral term */ 281c9681e5dSHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 282c9681e5dSHong Zhang } 283c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 284c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 285c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */ 286c9681e5dSHong Zhang if (ts->vec_costintegral) { 287c9681e5dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 288c9681e5dSHong Zhang } 289c9681e5dSHong Zhang } 290c9681e5dSHong Zhang 291c9681e5dSHong Zhang /* Build LHS for first-order adjoint */ 292c9681e5dSHong Zhang shift = 1./th->time_step; 293c9681e5dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 294c9681e5dSHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 295c9681e5dSHong Zhang 296c9681e5dSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 297c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 298c9681e5dSHong Zhang KSPConvergedReason kspreason; 299c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 300c9681e5dSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 301c9681e5dSHong Zhang if (kspreason < 0) { 302c9681e5dSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 303c9681e5dSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 304c9681e5dSHong Zhang } 305c9681e5dSHong Zhang } 306c9681e5dSHong Zhang 307c9681e5dSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 308c9681e5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 309c9681e5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 310c9681e5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 311c9681e5dSHong Zhang /* lambda_s^T F_UU w_1 */ 312c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 313c9681e5dSHong Zhang if (ts->vecs_fup) { 314c9681e5dSHong Zhang /* lambda_s^T F_UP w_2 */ 315c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 316c9681e5dSHong Zhang } 317c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 318c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 319c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 320c9681e5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 321c9681e5dSHong Zhang if (ts->vecs_fup) { 322c9681e5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 323c9681e5dSHong Zhang } 324c9681e5dSHong Zhang } 325c9681e5dSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 326c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 327*05c9950eSHong Zhang KSPConvergedReason kspreason; 328c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 329*05c9950eSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 330*05c9950eSHong Zhang if (kspreason < 0) { 331*05c9950eSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 332*05c9950eSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 333*05c9950eSHong Zhang } 334c9681e5dSHong Zhang } 335c9681e5dSHong Zhang } 336c9681e5dSHong Zhang 337c9681e5dSHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 338c9681e5dSHong Zhang shift = 0.0; 339c9681e5dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 340c9681e5dSHong Zhang ierr = MatScale(J,-1.);CHKERRQ(ierr); 341c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 342c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 343c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr); 344c9681e5dSHong Zhang ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 345c9681e5dSHong Zhang if (ts->vecs_sensi2) { 346c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 347c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr); 348c9681e5dSHong Zhang ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 349c9681e5dSHong Zhang } 350c9681e5dSHong Zhang } 351c9681e5dSHong Zhang if (ts->vecs_sensip) { 352c9681e5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */ 353c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 354c9681e5dSHong Zhang if (ts->vecs_fpu) { 355c9681e5dSHong Zhang /* lambda_s^T F_PU w_1 */ 356c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 357c9681e5dSHong Zhang } 358c9681e5dSHong Zhang if (ts->vecs_fpp) { 359c9681e5dSHong Zhang /* lambda_s^T F_PU w_2 */ 360c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 361c9681e5dSHong Zhang } 362c9681e5dSHong Zhang } 363c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 364c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 365c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 366c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 367c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 368c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 369c9681e5dSHong Zhang if (ts->vecs_fpu) { 370c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 371c9681e5dSHong Zhang } 372c9681e5dSHong Zhang if (ts->vecs_fpp) { 373c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 374c9681e5dSHong Zhang } 375c9681e5dSHong Zhang } 376c9681e5dSHong Zhang } 377c9681e5dSHong Zhang if (ts->vec_costintegral) { 378c9681e5dSHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 379c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 380c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 381c9681e5dSHong Zhang } 382c9681e5dSHong Zhang } 383c9681e5dSHong Zhang } 384c9681e5dSHong Zhang 385c9681e5dSHong Zhang if (ts->vecs_sensi2) { 386c9681e5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 387c9681e5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 388c9681e5dSHong Zhang } 389c9681e5dSHong Zhang th->status = TS_STEP_COMPLETE; 390c9681e5dSHong Zhang PetscFunctionReturn(0); 391c9681e5dSHong Zhang } 392c9681e5dSHong Zhang 39342f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 3942ca6e920SHong Zhang { 3952ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 396b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 397b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 3982ca6e920SHong Zhang PetscInt nadj; 3992ca6e920SHong Zhang Mat J,Jp; 4002ca6e920SHong Zhang KSP ksp; 4012ca6e920SHong Zhang PetscReal shift; 402b5b4017aSHong Zhang PetscScalar *xarr; 403b5b4017aSHong Zhang PetscErrorCode ierr; 4042ca6e920SHong Zhang 4052ca6e920SHong Zhang PetscFunctionBegin; 406c9681e5dSHong Zhang if (th->Theta == 1.) { 407c9681e5dSHong Zhang ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr); 408c9681e5dSHong Zhang PetscFunctionReturn(0); 409c9681e5dSHong Zhang } 4102ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 411302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 4122ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 413b5e0532dSHong Zhang 414b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 415022c081aSHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 416b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 417022c081aSHong Zhang th->time_step = -ts->time_step; 4182ca6e920SHong Zhang 419b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 4202c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 42105755b9cSHong Zhang if (th->endpoint) { 422c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 42305755b9cSHong Zhang } else { 424c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 42505755b9cSHong Zhang } 42605755b9cSHong Zhang } 427abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 428b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 429022c081aSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr); 4302c39e106SBarry Smith if (ts->vec_costintegral) { 431c9ad9de0SHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 43236eaed60SHong Zhang } 4332ca6e920SHong Zhang } 4343c54f92cSHong Zhang 435b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 436022c081aSHong Zhang shift = 1./(th->Theta*th->time_step); 4373c54f92cSHong Zhang if (th->endpoint) { 438022c081aSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 4393c54f92cSHong Zhang } else { 4403c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 4413c54f92cSHong Zhang } 4422ca6e920SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 4432ca6e920SHong Zhang 444b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 445abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 446b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 4472ca6e920SHong Zhang } 4483c54f92cSHong Zhang 449b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 450b5b4017aSHong Zhang if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 451b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 452b5b4017aSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 453b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 454b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 455b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 456b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 457b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 458b5b4017aSHong Zhang if (ts->vecs_fup) { 459b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 460b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 461b5b4017aSHong Zhang } 462b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 463b5b4017aSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 464b5b4017aSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],1./shift);CHKERRQ(ierr); 465b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 466b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 467b5b4017aSHong Zhang if (ts->vecs_fup) { 468b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 469b5b4017aSHong Zhang } 470b5b4017aSHong Zhang if (ts->vec_costintegral) { 471b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 472b5b4017aSHong Zhang } 473b5b4017aSHong Zhang } 474b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 475b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 476b5b4017aSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 477b5b4017aSHong Zhang } 478b5b4017aSHong Zhang } 479b5b4017aSHong Zhang 48036eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 4816affb6f8SHong Zhang if(th->endpoint) { /* two-stage Theta methods */ 482b5b4017aSHong Zhang if (th->Theta!=1.) { /* general case */ 483022c081aSHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 484b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 485b5b4017aSHong Zhang if (ts->vec_costintegral) { /* R_U at t_n */ 486c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 4878263dbbfSHong Zhang } 488abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 489b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 4903c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 491b5b4017aSHong Zhang if (ts->vec_costintegral) { 492b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1./shift,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 493b5b4017aSHong Zhang } 494b5b4017aSHong Zhang } 495b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* second-order */ 496b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 497b5b4017aSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 498b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 499b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 500b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 501b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 502b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 503b5b4017aSHong Zhang if (ts->vecs_fup) { 504b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 505b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 506b5b4017aSHong Zhang } 507b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 508b5b4017aSHong Zhang /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) R_U */ 509b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 510b5b4017aSHong Zhang ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr); 511b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 512b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 513b5b4017aSHong Zhang if (ts->vecs_fup) { 514b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fup[nadj]);CHKERRQ(ierr); 515b5b4017aSHong Zhang } 516b5b4017aSHong Zhang if (ts->vec_costintegral) { 517b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 518b5b4017aSHong Zhang } 519b5b4017aSHong Zhang } 5202ca6e920SHong Zhang } 521b5e0532dSHong Zhang } else { /* backward Euler */ 522b5e0532dSHong Zhang shift = 0.0; 523b5b4017aSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_u */ 524b5e0532dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 525b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 526022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 527b5b4017aSHong Zhang if (ts->vec_costintegral) { /* wrong? */ 528c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 529b5e0532dSHong Zhang } 530b5e0532dSHong Zhang } 531b5b4017aSHong Zhang if (ts->vecs_sensi2) { 532b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 533b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 534b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-th->time_step,VecsSensi2Temp[nadj]);CHKERRQ(ierr); 535b5b4017aSHong Zhang } 536b5b4017aSHong Zhang } 537b5e0532dSHong Zhang } 5383fd52205SHong Zhang 5393fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 540b5b4017aSHong Zhang /* U_{n+1} */ 541c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 542b5b4017aSHong Zhang if (ts->vec_costintegral) { 543b5b4017aSHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 544b5b4017aSHong Zhang } 545abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 546b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 547022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 5483fd52205SHong Zhang } 549b5b4017aSHong Zhang if (ts->vecs_sensip2) { /* second-order */ 550b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 551b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 552b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 553b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 554b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 555b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 556b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 557b5b4017aSHong Zhang if (ts->vecs_fpu) { 558b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 559b5b4017aSHong Zhang } 560b5b4017aSHong Zhang if (ts->vecs_fpp) { 561b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 562b5b4017aSHong Zhang } 563b5b4017aSHong Zhang if (ts->vec_costintegral) { 564b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 565b5b4017aSHong Zhang } 566b5b4017aSHong Zhang } 567b5b4017aSHong Zhang } 568b5b4017aSHong Zhang 569b5b4017aSHong Zhang /* U_s */ 570b5e0532dSHong Zhang if (th->Theta!=1.) { 571c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 572b5b4017aSHong Zhang if (ts->vec_costintegral) { 573b5b4017aSHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 574b5b4017aSHong Zhang } 575abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 576b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 577022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 578b5b4017aSHong Zhang if (ts->vecs_sensip2) { /* second-order */ 579b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 580b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 581b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 582b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 583b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 584b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 585b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 586b5b4017aSHong Zhang if (ts->vecs_fpu) { 587b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 588b5e0532dSHong Zhang } 589b5b4017aSHong Zhang if (ts->vecs_fpp) { 590b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 5913fd52205SHong Zhang } 5922c39e106SBarry Smith if (ts->vec_costintegral) { 593b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 59436eaed60SHong Zhang } 595b5b4017aSHong Zhang } 59636eaed60SHong Zhang } 59736eaed60SHong Zhang } 5983fd52205SHong Zhang } 599b5e0532dSHong Zhang } 6003fd52205SHong Zhang } else { /* one-stage case */ 6013c54f92cSHong Zhang shift = 0.0; 602a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 6032c39e106SBarry Smith if (ts->vec_costintegral) { 604c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 6053c54f92cSHong Zhang } 606abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 607b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 608022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 6092c39e106SBarry Smith if (ts->vec_costintegral) { 610c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 61136eaed60SHong Zhang } 6122ca6e920SHong Zhang } 6133fd52205SHong Zhang if (ts->vecs_sensip) { 614c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 615abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 616b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 617022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 6183fd52205SHong Zhang } 6192c39e106SBarry Smith if (ts->vec_costintegral) { 620c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 621abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 622022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 62336eaed60SHong Zhang } 62436eaed60SHong Zhang } 6253fd52205SHong Zhang } 6262ca6e920SHong Zhang } 6272ca6e920SHong Zhang 6282ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6292ca6e920SHong Zhang PetscFunctionReturn(0); 6302ca6e920SHong Zhang } 6312ca6e920SHong Zhang 632cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 633cd652676SJed Brown { 634cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 635be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 636cd652676SJed Brown PetscErrorCode ierr; 637cd652676SJed Brown 638cd652676SJed Brown PetscFunctionBegin; 639a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 640be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 641be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 642cd652676SJed Brown PetscFunctionReturn(0); 643cd652676SJed Brown } 644cd652676SJed Brown 6451566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 6461566a47fSLisandro Dalcin { 6471566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 6481566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 6491566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 6507453f775SEmil Constantinescu PetscReal wltea,wlter; 6511566a47fSLisandro Dalcin PetscErrorCode ierr; 6521566a47fSLisandro Dalcin 6531566a47fSLisandro Dalcin PetscFunctionBegin; 6542ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 6551566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 656fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 6571566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 658fecfb714SLisandro Dalcin { 659be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 660be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 6611566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 6621566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 6631566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 6641566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 6651566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 6667453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 6671566a47fSLisandro Dalcin } 6681566a47fSLisandro Dalcin if (order) *order = 2; 6691566a47fSLisandro Dalcin PetscFunctionReturn(0); 6701566a47fSLisandro Dalcin } 6711566a47fSLisandro Dalcin 6721566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 6731566a47fSLisandro Dalcin { 6741566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 67513b1b0ffSHong Zhang PetscInt ncost; 6761566a47fSLisandro Dalcin PetscErrorCode ierr; 6771566a47fSLisandro Dalcin 6781566a47fSLisandro Dalcin PetscFunctionBegin; 6791566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 6801566a47fSLisandro Dalcin if (ts->vec_costintegral && ts->costintegralfwd) { 6811566a47fSLisandro Dalcin ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 6821566a47fSLisandro Dalcin } 683715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 68413b1b0ffSHong Zhang if (ts->mat_sensip) { 68513b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 6866f92202bSHong Zhang } 6876f92202bSHong Zhang if (ts->vecs_integral_sensip) { 6886f92202bSHong Zhang for (ncost=0;ncost<ts->numcost;ncost++) { 6896f92202bSHong Zhang ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 6906f92202bSHong Zhang } 6916f92202bSHong Zhang } 692715f1b00SHong Zhang PetscFunctionReturn(0); 693715f1b00SHong Zhang } 694715f1b00SHong Zhang 695715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 696715f1b00SHong Zhang { 697715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 69813b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 699b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 70013b1b0ffSHong Zhang PetscInt ncost,ntlm; 701715f1b00SHong Zhang KSP ksp; 702715f1b00SHong Zhang Mat J,Jp; 703715f1b00SHong Zhang PetscReal shift; 70413b1b0ffSHong Zhang PetscScalar *barr,*xarr; 705715f1b00SHong Zhang PetscErrorCode ierr; 706715f1b00SHong Zhang 707715f1b00SHong Zhang PetscFunctionBegin; 70813b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 70913b1b0ffSHong Zhang 7106f92202bSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 7116f92202bSHong Zhang if (ts->vecs_integral_sensip) { 7126f92202bSHong Zhang ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 7136f92202bSHong Zhang } 7146f92202bSHong Zhang } 7156f92202bSHong Zhang 716715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 717715f1b00SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 718715f1b00SHong Zhang 719715f1b00SHong Zhang /* Build RHS */ 720715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 721715f1b00SHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 722715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 72313b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 72413b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 725715f1b00SHong Zhang 726022c081aSHong Zhang /* Add the f_p forcing terms */ 7270e11c21fSHong Zhang if (ts->Jacp) { 728c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 72913b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(1.-th->Theta)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 730c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 73113b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 7320e11c21fSHong Zhang } 733715f1b00SHong Zhang } else { /* 1-stage method */ 734715f1b00SHong Zhang shift = 0.0; 735715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 73613b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 73713b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 73813b1b0ffSHong Zhang 739022c081aSHong Zhang /* Add the f_p forcing terms */ 7400e11c21fSHong Zhang if (ts->Jacp) { 741c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 74213b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 743715f1b00SHong Zhang } 7440e11c21fSHong Zhang } 745715f1b00SHong Zhang 746715f1b00SHong Zhang /* Build LHS */ 747715f1b00SHong Zhang shift = 1/(th->Theta*th->time_step); 748715f1b00SHong Zhang if (th->endpoint) { 749715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 750715f1b00SHong Zhang } else { 751715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 752715f1b00SHong Zhang } 753715f1b00SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 754715f1b00SHong Zhang 755715f1b00SHong Zhang /* 756715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 757c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 758715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 759715f1b00SHong Zhang */ 760715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 761715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 762c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 7630e11c21fSHong Zhang if (ts->vecs_drdp) { 764c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 7650e11c21fSHong Zhang } 766715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 767c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 7680e11c21fSHong Zhang if (ts->vecs_drdp) { 76913b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 7700e11c21fSHong Zhang } 771715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 772715f1b00SHong Zhang } 773715f1b00SHong Zhang } 774715f1b00SHong Zhang } 775715f1b00SHong Zhang 776715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 777022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 778b5b4017aSHong Zhang KSPConvergedReason kspreason; 77913b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 780b5b4017aSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 781715f1b00SHong Zhang if (th->endpoint) { 78213b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 783b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 784b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 785b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 78613b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 787715f1b00SHong Zhang } else { 788b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 789715f1b00SHong Zhang } 790b5b4017aSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 791b5b4017aSHong Zhang if (kspreason < 0) { 792b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 793b5b4017aSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 794b5b4017aSHong Zhang } 795b5b4017aSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 79613b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 797715f1b00SHong Zhang } 798715f1b00SHong Zhang 799b5b4017aSHong Zhang 80013b1b0ffSHong Zhang /* 80113b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 802c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 80313b1b0ffSHong Zhang */ 80413b1b0ffSHong Zhang if (ts->vecs_integral_sensip) { 80513b1b0ffSHong Zhang if (!th->endpoint) { 80613b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 807c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 8080e11c21fSHong Zhang if (ts->vecs_drdp) { 809c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 8100e11c21fSHong Zhang } 81113b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 812c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 8130e11c21fSHong Zhang if (ts->vecs_drdp) { 81413b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 8150e11c21fSHong Zhang } 81613b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 81713b1b0ffSHong Zhang } 81813b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 81913b1b0ffSHong Zhang } else { 820c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 8210e11c21fSHong Zhang if (ts->vecs_drdp) { 822c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 8230e11c21fSHong Zhang } 82413b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 825c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 8260e11c21fSHong Zhang if (ts->vecs_drdp) { 82713b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 8280e11c21fSHong Zhang } 82913b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 83013b1b0ffSHong Zhang } 83113b1b0ffSHong Zhang } 83213b1b0ffSHong Zhang } else { 83313b1b0ffSHong Zhang if (!th->endpoint) { 83413b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 835715f1b00SHong Zhang } 836715f1b00SHong Zhang } 8371566a47fSLisandro Dalcin PetscFunctionReturn(0); 8381566a47fSLisandro Dalcin } 8391566a47fSLisandro Dalcin 8406affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 8416affb6f8SHong Zhang { 8426affb6f8SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 8436affb6f8SHong Zhang 8446affb6f8SHong Zhang PetscFunctionBegin; 8456affb6f8SHong Zhang if (ns) *ns = 1; 8466affb6f8SHong Zhang if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 8476affb6f8SHong Zhang PetscFunctionReturn(0); 8486affb6f8SHong Zhang } 8496affb6f8SHong Zhang 850316643e7SJed Brown /*------------------------------------------------------------*/ 851277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 852316643e7SJed Brown { 853316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 854316643e7SJed Brown PetscErrorCode ierr; 855316643e7SJed Brown 856316643e7SJed Brown PetscFunctionBegin; 8576bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 8586bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 8593b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 860eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 8611566a47fSLisandro Dalcin 8621566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 8631566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 8641566a47fSLisandro Dalcin 865b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 866715f1b00SHong Zhang if (ts->forward_solve) { 867715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 868715f1b00SHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 8696f92202bSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 870715f1b00SHong Zhang } 871b5b4017aSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 87213b1b0ffSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 87313b1b0ffSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 874715f1b00SHong Zhang } 875b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 876b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 877b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 878b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 879b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 880b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 881715f1b00SHong Zhang 882277b19d0SLisandro Dalcin PetscFunctionReturn(0); 883277b19d0SLisandro Dalcin } 884277b19d0SLisandro Dalcin 885ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts) 886ece46509SHong Zhang { 887ece46509SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 888ece46509SHong Zhang PetscErrorCode ierr; 889ece46509SHong Zhang 890ece46509SHong Zhang PetscFunctionBegin; 891ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 892ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 893ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 894ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 895ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 896ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 897ece46509SHong Zhang PetscFunctionReturn(0); 898ece46509SHong Zhang } 899ece46509SHong Zhang 900277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 901277b19d0SLisandro Dalcin { 902277b19d0SLisandro Dalcin PetscErrorCode ierr; 903277b19d0SLisandro Dalcin 904277b19d0SLisandro Dalcin PetscFunctionBegin; 905277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 906b3a6b972SJed Brown if (ts->dm) { 907b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 908b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 909b3a6b972SJed Brown } 910277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 911bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 912bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 913bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 914bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 915316643e7SJed Brown PetscFunctionReturn(0); 916316643e7SJed Brown } 917316643e7SJed Brown 918316643e7SJed Brown /* 919316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9202b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 921316643e7SJed Brown */ 9220f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 923316643e7SJed Brown { 924316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 925316643e7SJed Brown PetscErrorCode ierr; 9267445fe48SJed Brown Vec X0,Xdot; 9277445fe48SJed Brown DM dm,dmsave; 9281566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 929316643e7SJed Brown 930316643e7SJed Brown PetscFunctionBegin; 9317445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 9325a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9337445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 934b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 9357445fe48SJed Brown 9367445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9377445fe48SJed Brown dmsave = ts->dm; 9387445fe48SJed Brown ts->dm = dm; 9397445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 9407445fe48SJed Brown ts->dm = dmsave; 9410d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 942316643e7SJed Brown PetscFunctionReturn(0); 943316643e7SJed Brown } 944316643e7SJed Brown 945d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 946316643e7SJed Brown { 947316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 948316643e7SJed Brown PetscErrorCode ierr; 9497445fe48SJed Brown Vec Xdot; 9507445fe48SJed Brown DM dm,dmsave; 9511566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 952316643e7SJed Brown 953316643e7SJed Brown PetscFunctionBegin; 9547445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 955be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 9560298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 9577445fe48SJed Brown 9587445fe48SJed Brown dmsave = ts->dm; 9597445fe48SJed Brown ts->dm = dm; 960d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 9617445fe48SJed Brown ts->dm = dmsave; 9620298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 963316643e7SJed Brown PetscFunctionReturn(0); 964316643e7SJed Brown } 965316643e7SJed Brown 966715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 967715f1b00SHong Zhang { 968715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 969715f1b00SHong Zhang PetscErrorCode ierr; 970715f1b00SHong Zhang 971715f1b00SHong Zhang PetscFunctionBegin; 972022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 97313b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 97413b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 975715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 976715f1b00SHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 977715f1b00SHong Zhang } 9786f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 97913b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 98013b1b0ffSHong Zhang 9816f92202bSHong Zhang if (ts->vecs_integral_sensip) { 9826f92202bSHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 9836f92202bSHong Zhang } 984b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 985b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 986715f1b00SHong Zhang PetscFunctionReturn(0); 987715f1b00SHong Zhang } 988715f1b00SHong Zhang 9897adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts) 9907adebddeSHong Zhang { 9917adebddeSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 9927adebddeSHong Zhang PetscErrorCode ierr; 9937adebddeSHong Zhang 9947adebddeSHong Zhang PetscFunctionBegin; 9957adebddeSHong Zhang if (ts->vecs_integral_sensip) { 9967adebddeSHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 9977adebddeSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 9987adebddeSHong Zhang } 9997adebddeSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 10007adebddeSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 10017adebddeSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 10027adebddeSHong Zhang PetscFunctionReturn(0); 10037adebddeSHong Zhang } 10047adebddeSHong Zhang 1005316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 1006316643e7SJed Brown { 1007316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 10082ffb9264SLisandro Dalcin PetscBool match; 1009316643e7SJed Brown PetscErrorCode ierr; 1010316643e7SJed Brown 1011316643e7SJed Brown PetscFunctionBegin; 1012b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 1013b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 1014b1cb36f3SHong Zhang } 101539d13666SHong Zhang if (!th->X) { 1016316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 101739d13666SHong Zhang } 101839d13666SHong Zhang if (!th->Xdot) { 1019316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 102039d13666SHong Zhang } 102139d13666SHong Zhang if (!th->X0) { 10223b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 102339d13666SHong Zhang } 10241566a47fSLisandro Dalcin if (th->endpoint) { 10251566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 10267445fe48SJed Brown } 10273b1890cdSShri Abhyankar 10281566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 10291566a47fSLisandro Dalcin 10301566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 10311566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10321566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10331566a47fSLisandro Dalcin 10341566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 10351566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 10362ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 10372ffb9264SLisandro Dalcin if (!match) { 10381566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 10391566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 10403b1890cdSShri Abhyankar } 10411566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 1042b4dd3bf9SBarry Smith PetscFunctionReturn(0); 1043b4dd3bf9SBarry Smith } 10440c7376e5SHong Zhang 1045b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1046b4dd3bf9SBarry Smith 1047b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1048b4dd3bf9SBarry Smith { 1049b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 1050b4dd3bf9SBarry Smith PetscErrorCode ierr; 1051b4dd3bf9SBarry Smith 1052b4dd3bf9SBarry Smith PetscFunctionBegin; 1053b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 1054b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 10552ca6e920SHong Zhang if (ts->vecs_sensip) { 1056b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 10572ca6e920SHong Zhang } 1058b5b4017aSHong Zhang if (ts->vecs_sensi2) { 1059b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 1060b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 1061b5b4017aSHong Zhang } 1062c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 1063c9681e5dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 1064b5b4017aSHong Zhang } 1065316643e7SJed Brown PetscFunctionReturn(0); 1066316643e7SJed Brown } 1067316643e7SJed Brown 10684416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 1069316643e7SJed Brown { 1070316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1071316643e7SJed Brown PetscErrorCode ierr; 1072316643e7SJed Brown 1073316643e7SJed Brown PetscFunctionBegin; 1074e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 1075316643e7SJed Brown { 10760298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 10770298fd71SBarry 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); 107803842d09SLisandro 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); 1079316643e7SJed Brown } 1080316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1081316643e7SJed Brown PetscFunctionReturn(0); 1082316643e7SJed Brown } 1083316643e7SJed Brown 1084316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 1085316643e7SJed Brown { 1086316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1087ace3abfcSBarry Smith PetscBool iascii; 1088316643e7SJed Brown PetscErrorCode ierr; 1089316643e7SJed Brown 1090316643e7SJed Brown PetscFunctionBegin; 1091251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1092316643e7SJed Brown if (iascii) { 10937c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 1094ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1095316643e7SJed Brown } 1096316643e7SJed Brown PetscFunctionReturn(0); 1097316643e7SJed Brown } 1098316643e7SJed Brown 1099560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 11000de4c49aSJed Brown { 11010de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11020de4c49aSJed Brown 11030de4c49aSJed Brown PetscFunctionBegin; 11040de4c49aSJed Brown *theta = th->Theta; 11050de4c49aSJed Brown PetscFunctionReturn(0); 11060de4c49aSJed Brown } 11070de4c49aSJed Brown 1108560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 11090de4c49aSJed Brown { 11100de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11110de4c49aSJed Brown 11120de4c49aSJed Brown PetscFunctionBegin; 11137c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 11140de4c49aSJed Brown th->Theta = theta; 11151566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 11160de4c49aSJed Brown PetscFunctionReturn(0); 11170de4c49aSJed Brown } 1118eb284becSJed Brown 1119560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 112026f2ff8fSLisandro Dalcin { 112126f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 112226f2ff8fSLisandro Dalcin 112326f2ff8fSLisandro Dalcin PetscFunctionBegin; 112426f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 112526f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 112626f2ff8fSLisandro Dalcin } 112726f2ff8fSLisandro Dalcin 1128560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 1129eb284becSJed Brown { 1130eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1131eb284becSJed Brown 1132eb284becSJed Brown PetscFunctionBegin; 1133eb284becSJed Brown th->endpoint = flg; 1134eb284becSJed Brown PetscFunctionReturn(0); 1135eb284becSJed Brown } 11360de4c49aSJed Brown 1137f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1138f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 1139f9c1d6abSBarry Smith { 1140f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 1141f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 11423fd8ae06SJed Brown const PetscReal one = 1.0; 1143f9c1d6abSBarry Smith 1144f9c1d6abSBarry Smith PetscFunctionBegin; 11453fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 1146f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 1147f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 1148f9c1d6abSBarry Smith PetscFunctionReturn(0); 1149f9c1d6abSBarry Smith } 1150f9c1d6abSBarry Smith #endif 1151f9c1d6abSBarry Smith 115242682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 115342682096SHong Zhang { 115442682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 115542682096SHong Zhang 115642682096SHong Zhang PetscFunctionBegin; 11571566a47fSLisandro Dalcin if (ns) *ns = 1; 11581566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 115942682096SHong Zhang PetscFunctionReturn(0); 116042682096SHong Zhang } 1161f9c1d6abSBarry Smith 1162316643e7SJed Brown /* ------------------------------------------------------------ */ 1163316643e7SJed Brown /*MC 116496f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1165316643e7SJed Brown 1166316643e7SJed Brown Level: beginner 1167316643e7SJed Brown 11684eb428fdSBarry Smith Options Database: 1169c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1170c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 117103842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 11724eb428fdSBarry Smith 1173eb284becSJed Brown Notes: 11740c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 11750c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 11764eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 11774eb428fdSBarry Smith 1178eb284becSJed Brown This method can be applied to DAE. 1179eb284becSJed Brown 1180eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 1181eb284becSJed Brown 1182eb284becSJed Brown .vb 1183eb284becSJed Brown Theta | Theta 1184eb284becSJed Brown ------------- 1185eb284becSJed Brown | 1 1186eb284becSJed Brown .ve 1187eb284becSJed Brown 1188eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1189eb284becSJed Brown 1190eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1191eb284becSJed Brown 1192eb284becSJed Brown .vb 1193eb284becSJed Brown 0 | 0 0 1194eb284becSJed Brown 1 | 1-Theta Theta 1195eb284becSJed Brown ------------------- 1196eb284becSJed Brown | 1-Theta Theta 1197eb284becSJed Brown .ve 1198eb284becSJed Brown 1199eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1200eb284becSJed Brown 1201eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1202eb284becSJed Brown 1203eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1204eb284becSJed Brown 12054eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1206eb284becSJed Brown 1207eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1208316643e7SJed Brown 1209316643e7SJed Brown M*/ 12108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1211316643e7SJed Brown { 1212316643e7SJed Brown TS_Theta *th; 1213316643e7SJed Brown PetscErrorCode ierr; 1214316643e7SJed Brown 1215316643e7SJed Brown PetscFunctionBegin; 1216277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1217ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1218316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1219316643e7SJed Brown ts->ops->view = TSView_Theta; 1220316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 122142f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1222ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1223316643e7SJed Brown ts->ops->step = TSStep_Theta; 1224cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12251566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 122624655328SShri ts->ops->rollback = TSRollBack_Theta; 1227316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12280f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12290f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1230f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1231f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1232f9c1d6abSBarry Smith #endif 123342682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 123442f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1235b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1236b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12372ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12386affb6f8SHong Zhang 1239715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12407adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1241715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12426affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1243316643e7SJed Brown 1244efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1245efd4aadfSBarry Smith 1246b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1247316643e7SJed Brown ts->data = (void*)th; 1248316643e7SJed Brown 1249715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1250715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1251715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1252b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1253715f1b00SHong Zhang 12546f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1255316643e7SJed Brown th->Theta = 0.5; 12561566a47fSLisandro Dalcin th->order = 2; 1257bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1258bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1259bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1260bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1261316643e7SJed Brown PetscFunctionReturn(0); 1262316643e7SJed Brown } 12630de4c49aSJed Brown 12640de4c49aSJed Brown /*@ 12650de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 12660de4c49aSJed Brown 12670de4c49aSJed Brown Not Collective 12680de4c49aSJed Brown 12690de4c49aSJed Brown Input Parameter: 12700de4c49aSJed Brown . ts - timestepping context 12710de4c49aSJed Brown 12720de4c49aSJed Brown Output Parameter: 12730de4c49aSJed Brown . theta - stage abscissa 12740de4c49aSJed Brown 12750de4c49aSJed Brown Note: 12760de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 12770de4c49aSJed Brown 12780de4c49aSJed Brown Level: Advanced 12790de4c49aSJed Brown 12800de4c49aSJed Brown .seealso: TSThetaSetTheta() 12810de4c49aSJed Brown @*/ 12827087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 12830de4c49aSJed Brown { 12844ac538c5SBarry Smith PetscErrorCode ierr; 12850de4c49aSJed Brown 12860de4c49aSJed Brown PetscFunctionBegin; 1287afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12880de4c49aSJed Brown PetscValidPointer(theta,2); 12894ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 12900de4c49aSJed Brown PetscFunctionReturn(0); 12910de4c49aSJed Brown } 12920de4c49aSJed Brown 12930de4c49aSJed Brown /*@ 12940de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 12950de4c49aSJed Brown 12960de4c49aSJed Brown Not Collective 12970de4c49aSJed Brown 12980de4c49aSJed Brown Input Parameter: 12990de4c49aSJed Brown + ts - timestepping context 13000de4c49aSJed Brown - theta - stage abscissa 13010de4c49aSJed Brown 13020de4c49aSJed Brown Options Database: 13030de4c49aSJed Brown . -ts_theta_theta <theta> 13040de4c49aSJed Brown 13050de4c49aSJed Brown Level: Intermediate 13060de4c49aSJed Brown 13070de4c49aSJed Brown .seealso: TSThetaGetTheta() 13080de4c49aSJed Brown @*/ 13097087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 13100de4c49aSJed Brown { 13114ac538c5SBarry Smith PetscErrorCode ierr; 13120de4c49aSJed Brown 13130de4c49aSJed Brown PetscFunctionBegin; 1314afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13154ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 13160de4c49aSJed Brown PetscFunctionReturn(0); 13170de4c49aSJed Brown } 1318f33bbcb6SJed Brown 131926f2ff8fSLisandro Dalcin /*@ 132026f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 132126f2ff8fSLisandro Dalcin 132226f2ff8fSLisandro Dalcin Not Collective 132326f2ff8fSLisandro Dalcin 132426f2ff8fSLisandro Dalcin Input Parameter: 132526f2ff8fSLisandro Dalcin . ts - timestepping context 132626f2ff8fSLisandro Dalcin 132726f2ff8fSLisandro Dalcin Output Parameter: 132826f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 132926f2ff8fSLisandro Dalcin 133026f2ff8fSLisandro Dalcin Level: Advanced 133126f2ff8fSLisandro Dalcin 133226f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 133326f2ff8fSLisandro Dalcin @*/ 133426f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 133526f2ff8fSLisandro Dalcin { 133626f2ff8fSLisandro Dalcin PetscErrorCode ierr; 133726f2ff8fSLisandro Dalcin 133826f2ff8fSLisandro Dalcin PetscFunctionBegin; 133926f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 134026f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1341163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 134226f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 134326f2ff8fSLisandro Dalcin } 134426f2ff8fSLisandro Dalcin 1345eb284becSJed Brown /*@ 1346eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1347eb284becSJed Brown 1348eb284becSJed Brown Not Collective 1349eb284becSJed Brown 1350eb284becSJed Brown Input Parameter: 1351eb284becSJed Brown + ts - timestepping context 1352eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1353eb284becSJed Brown 1354eb284becSJed Brown Options Database: 1355eb284becSJed Brown . -ts_theta_endpoint <flg> 1356eb284becSJed Brown 1357eb284becSJed Brown Level: Intermediate 1358eb284becSJed Brown 1359eb284becSJed Brown .seealso: TSTHETA, TSCN 1360eb284becSJed Brown @*/ 1361eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1362eb284becSJed Brown { 1363eb284becSJed Brown PetscErrorCode ierr; 1364eb284becSJed Brown 1365eb284becSJed Brown PetscFunctionBegin; 1366eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1367eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1368eb284becSJed Brown PetscFunctionReturn(0); 1369eb284becSJed Brown } 1370eb284becSJed Brown 1371f33bbcb6SJed Brown /* 1372f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1373f33bbcb6SJed Brown * The creation functions for these specializations are below. 1374f33bbcb6SJed Brown */ 1375f33bbcb6SJed Brown 13761566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 13771566a47fSLisandro Dalcin { 13781566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 13791566a47fSLisandro Dalcin PetscErrorCode ierr; 13801566a47fSLisandro Dalcin 13811566a47fSLisandro Dalcin PetscFunctionBegin; 13821566a47fSLisandro 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"); 13831566a47fSLisandro 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"); 13841566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 13851566a47fSLisandro Dalcin PetscFunctionReturn(0); 13861566a47fSLisandro Dalcin } 13871566a47fSLisandro Dalcin 1388f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1389f33bbcb6SJed Brown { 1390f33bbcb6SJed Brown PetscFunctionBegin; 1391f33bbcb6SJed Brown PetscFunctionReturn(0); 1392f33bbcb6SJed Brown } 1393f33bbcb6SJed Brown 1394f33bbcb6SJed Brown /*MC 1395f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1396f33bbcb6SJed Brown 1397f33bbcb6SJed Brown Level: beginner 1398f33bbcb6SJed Brown 13994eb428fdSBarry Smith Notes: 1400c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 14014eb428fdSBarry Smith 14021566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 14034eb428fdSBarry Smith 1404f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1405f33bbcb6SJed Brown 1406f33bbcb6SJed Brown M*/ 14078cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1408f33bbcb6SJed Brown { 1409f33bbcb6SJed Brown PetscErrorCode ierr; 1410f33bbcb6SJed Brown 1411f33bbcb6SJed Brown PetscFunctionBegin; 1412f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1413f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 14141566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 14150c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1416f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1417f33bbcb6SJed Brown PetscFunctionReturn(0); 1418f33bbcb6SJed Brown } 1419f33bbcb6SJed Brown 14201566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 14211566a47fSLisandro Dalcin { 14221566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14231566a47fSLisandro Dalcin PetscErrorCode ierr; 14241566a47fSLisandro Dalcin 14251566a47fSLisandro Dalcin PetscFunctionBegin; 14261566a47fSLisandro 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"); 14271566a47fSLisandro 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"); 14281566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14291566a47fSLisandro Dalcin PetscFunctionReturn(0); 14301566a47fSLisandro Dalcin } 14311566a47fSLisandro Dalcin 1432f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1433f33bbcb6SJed Brown { 1434f33bbcb6SJed Brown PetscFunctionBegin; 1435f33bbcb6SJed Brown PetscFunctionReturn(0); 1436f33bbcb6SJed Brown } 1437f33bbcb6SJed Brown 1438f33bbcb6SJed Brown /*MC 1439f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1440f33bbcb6SJed Brown 1441f33bbcb6SJed Brown Level: beginner 1442f33bbcb6SJed Brown 1443f33bbcb6SJed Brown Notes: 14447cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 14457cf5af47SJed Brown 14467cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1447f33bbcb6SJed Brown 1448f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1449f33bbcb6SJed Brown 1450f33bbcb6SJed Brown M*/ 14518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1452f33bbcb6SJed Brown { 1453f33bbcb6SJed Brown PetscErrorCode ierr; 1454f33bbcb6SJed Brown 1455f33bbcb6SJed Brown PetscFunctionBegin; 1456f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1457f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1458eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 14590c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1460f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1461f33bbcb6SJed Brown PetscFunctionReturn(0); 1462f33bbcb6SJed Brown } 1463