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++) { 32705c9950eSHong Zhang KSPConvergedReason kspreason; 328c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 32905c9950eSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 33005c9950eSHong Zhang if (kspreason < 0) { 33105c9950eSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 33205c9950eSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 33305c9950eSHong 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++) { 446*1d14d03bSHong Zhang KSPConvergedReason kspreason; 447b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 448*1d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 449*1d14d03bSHong Zhang if (kspreason < 0) { 450*1d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 451*1d14d03bSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 452*1d14d03bSHong Zhang } 4532ca6e920SHong Zhang } 4543c54f92cSHong Zhang 455*1d14d03bSHong Zhang /* Second-order adjoint */ 456b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 457b5b4017aSHong Zhang if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 458b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 459b5b4017aSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 460b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 461b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 462b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 463b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 464b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 465b5b4017aSHong Zhang if (ts->vecs_fup) { 466b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 467b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 468b5b4017aSHong Zhang } 469b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 470b5b4017aSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 471*1d14d03bSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 472b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 473b5b4017aSHong Zhang if (ts->vecs_fup) { 474b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 475b5b4017aSHong Zhang } 476b5b4017aSHong Zhang } 477b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 478b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 479*1d14d03bSHong Zhang KSPConvergedReason kspreason; 480*1d14d03bSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 481*1d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 482*1d14d03bSHong Zhang if (kspreason < 0) { 483*1d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 484*1d14d03bSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 485*1d14d03bSHong Zhang } 486b5b4017aSHong Zhang } 487b5b4017aSHong Zhang } 488b5b4017aSHong Zhang 48936eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 490*1d14d03bSHong Zhang if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 491022c081aSHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 492b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 493b5b4017aSHong Zhang if (ts->vec_costintegral) { /* R_U at t_n */ 494c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 4958263dbbfSHong Zhang } 496abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 497b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 498b5b4017aSHong Zhang if (ts->vec_costintegral) { 499*1d14d03bSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 500b5b4017aSHong Zhang } 501*1d14d03bSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 502b5b4017aSHong Zhang } 503*1d14d03bSHong Zhang 504*1d14d03bSHong Zhang /* Second-order adjoint */ 505*1d14d03bSHong Zhang if (ts->vecs_sensi2) { /* U_n */ 506b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 507b5b4017aSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 508b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 509b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 510b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 511b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 512b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 513b5b4017aSHong Zhang if (ts->vecs_fup) { 514b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 515b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 516b5b4017aSHong Zhang } 517b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 518*1d14d03bSHong 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_UU w_2 */ 519b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 520*1d14d03bSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 521b5b4017aSHong Zhang if (ts->vecs_fup) { 522*1d14d03bSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 523b5b4017aSHong Zhang } 524*1d14d03bSHong Zhang ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr); 525b5b4017aSHong Zhang } 526b5e0532dSHong Zhang } 5273fd52205SHong Zhang 5283fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 529b5b4017aSHong Zhang /* U_{n+1} */ 530c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 531b5b4017aSHong Zhang if (ts->vec_costintegral) { 532b5b4017aSHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 533b5b4017aSHong Zhang } 534abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 535b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 536022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 5373fd52205SHong Zhang } 538b5b4017aSHong Zhang if (ts->vecs_sensip2) { /* second-order */ 539b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 540b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 541b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 542b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 543b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 544b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 545b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 546b5b4017aSHong Zhang if (ts->vecs_fpu) { 547*1d14d03bSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 548b5b4017aSHong Zhang } 549b5b4017aSHong Zhang if (ts->vecs_fpp) { 550*1d14d03bSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 551b5b4017aSHong Zhang } 552b5b4017aSHong Zhang } 553b5b4017aSHong Zhang } 554b5b4017aSHong Zhang 555b5b4017aSHong Zhang /* U_s */ 556c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 557b5b4017aSHong Zhang if (ts->vec_costintegral) { 558b5b4017aSHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 559b5b4017aSHong Zhang } 560abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 561b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 562022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 563b5b4017aSHong Zhang if (ts->vecs_sensip2) { /* second-order */ 564b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 565b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 566b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 567b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 568b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 569b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 570b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 571b5b4017aSHong Zhang if (ts->vecs_fpu) { 572*1d14d03bSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 573b5e0532dSHong Zhang } 574b5b4017aSHong Zhang if (ts->vecs_fpp) { 575*1d14d03bSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 576b5b4017aSHong Zhang } 57736eaed60SHong Zhang } 57836eaed60SHong Zhang } 5793fd52205SHong Zhang } 580b5e0532dSHong Zhang } 5813fd52205SHong Zhang } else { /* one-stage case */ 5823c54f92cSHong Zhang shift = 0.0; 583a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 5842c39e106SBarry Smith if (ts->vec_costintegral) { 585c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 5863c54f92cSHong Zhang } 587abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 588b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 589022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 5902c39e106SBarry Smith if (ts->vec_costintegral) { 591c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 59236eaed60SHong Zhang } 5932ca6e920SHong Zhang } 5943fd52205SHong Zhang if (ts->vecs_sensip) { 595c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 596abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 597b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 598022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 5993fd52205SHong Zhang } 6002c39e106SBarry Smith if (ts->vec_costintegral) { 601c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 602abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 603022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 60436eaed60SHong Zhang } 60536eaed60SHong Zhang } 6063fd52205SHong Zhang } 6072ca6e920SHong Zhang } 6082ca6e920SHong Zhang 6092ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6102ca6e920SHong Zhang PetscFunctionReturn(0); 6112ca6e920SHong Zhang } 6122ca6e920SHong Zhang 613cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 614cd652676SJed Brown { 615cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 616be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 617cd652676SJed Brown PetscErrorCode ierr; 618cd652676SJed Brown 619cd652676SJed Brown PetscFunctionBegin; 620a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 621be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 622be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 623cd652676SJed Brown PetscFunctionReturn(0); 624cd652676SJed Brown } 625cd652676SJed Brown 6261566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 6271566a47fSLisandro Dalcin { 6281566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 6291566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 6301566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 6317453f775SEmil Constantinescu PetscReal wltea,wlter; 6321566a47fSLisandro Dalcin PetscErrorCode ierr; 6331566a47fSLisandro Dalcin 6341566a47fSLisandro Dalcin PetscFunctionBegin; 6352ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 6361566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 637fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 6381566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 639fecfb714SLisandro Dalcin { 640be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 641be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 6421566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 6431566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 6441566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 6451566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 6461566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 6477453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 6481566a47fSLisandro Dalcin } 6491566a47fSLisandro Dalcin if (order) *order = 2; 6501566a47fSLisandro Dalcin PetscFunctionReturn(0); 6511566a47fSLisandro Dalcin } 6521566a47fSLisandro Dalcin 6531566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 6541566a47fSLisandro Dalcin { 6551566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 65613b1b0ffSHong Zhang PetscInt ncost; 6571566a47fSLisandro Dalcin PetscErrorCode ierr; 6581566a47fSLisandro Dalcin 6591566a47fSLisandro Dalcin PetscFunctionBegin; 6601566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 6611566a47fSLisandro Dalcin if (ts->vec_costintegral && ts->costintegralfwd) { 6621566a47fSLisandro Dalcin ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 6631566a47fSLisandro Dalcin } 664715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 66513b1b0ffSHong Zhang if (ts->mat_sensip) { 66613b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 6676f92202bSHong Zhang } 6686f92202bSHong Zhang if (ts->vecs_integral_sensip) { 6696f92202bSHong Zhang for (ncost=0;ncost<ts->numcost;ncost++) { 6706f92202bSHong Zhang ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 6716f92202bSHong Zhang } 6726f92202bSHong Zhang } 673715f1b00SHong Zhang PetscFunctionReturn(0); 674715f1b00SHong Zhang } 675715f1b00SHong Zhang 676715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 677715f1b00SHong Zhang { 678715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 67913b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 680b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 68113b1b0ffSHong Zhang PetscInt ncost,ntlm; 682715f1b00SHong Zhang KSP ksp; 683715f1b00SHong Zhang Mat J,Jp; 684715f1b00SHong Zhang PetscReal shift; 68513b1b0ffSHong Zhang PetscScalar *barr,*xarr; 686715f1b00SHong Zhang PetscErrorCode ierr; 687715f1b00SHong Zhang 688715f1b00SHong Zhang PetscFunctionBegin; 68913b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 69013b1b0ffSHong Zhang 6916f92202bSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 6926f92202bSHong Zhang if (ts->vecs_integral_sensip) { 6936f92202bSHong Zhang ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 6946f92202bSHong Zhang } 6956f92202bSHong Zhang } 6966f92202bSHong Zhang 697715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 698715f1b00SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 699715f1b00SHong Zhang 700715f1b00SHong Zhang /* Build RHS */ 701715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 702715f1b00SHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 703715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 70413b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 70513b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 706715f1b00SHong Zhang 707022c081aSHong Zhang /* Add the f_p forcing terms */ 7080e11c21fSHong Zhang if (ts->Jacp) { 709c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 71013b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(1.-th->Theta)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 711c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 71213b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 7130e11c21fSHong Zhang } 714715f1b00SHong Zhang } else { /* 1-stage method */ 715715f1b00SHong Zhang shift = 0.0; 716715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 71713b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 71813b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 71913b1b0ffSHong Zhang 720022c081aSHong Zhang /* Add the f_p forcing terms */ 7210e11c21fSHong Zhang if (ts->Jacp) { 722c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 72313b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 724715f1b00SHong Zhang } 7250e11c21fSHong Zhang } 726715f1b00SHong Zhang 727715f1b00SHong Zhang /* Build LHS */ 728715f1b00SHong Zhang shift = 1/(th->Theta*th->time_step); 729715f1b00SHong Zhang if (th->endpoint) { 730715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 731715f1b00SHong Zhang } else { 732715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 733715f1b00SHong Zhang } 734715f1b00SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 735715f1b00SHong Zhang 736715f1b00SHong Zhang /* 737715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 738c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 739715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 740715f1b00SHong Zhang */ 741715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 742715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 743c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 7440e11c21fSHong Zhang if (ts->vecs_drdp) { 745c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 7460e11c21fSHong Zhang } 747715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 748c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 7490e11c21fSHong Zhang if (ts->vecs_drdp) { 75013b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 7510e11c21fSHong Zhang } 752715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 753715f1b00SHong Zhang } 754715f1b00SHong Zhang } 755715f1b00SHong Zhang } 756715f1b00SHong Zhang 757715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 758022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 759b5b4017aSHong Zhang KSPConvergedReason kspreason; 76013b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 761b5b4017aSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 762715f1b00SHong Zhang if (th->endpoint) { 76313b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 764b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 765b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 766b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 76713b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 768715f1b00SHong Zhang } else { 769b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 770715f1b00SHong Zhang } 771b5b4017aSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 772b5b4017aSHong Zhang if (kspreason < 0) { 773b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 774b5b4017aSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 775b5b4017aSHong Zhang } 776b5b4017aSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 77713b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 778715f1b00SHong Zhang } 779715f1b00SHong Zhang 780b5b4017aSHong Zhang 78113b1b0ffSHong Zhang /* 78213b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 783c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 78413b1b0ffSHong Zhang */ 78513b1b0ffSHong Zhang if (ts->vecs_integral_sensip) { 78613b1b0ffSHong Zhang if (!th->endpoint) { 78713b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 788c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 7890e11c21fSHong Zhang if (ts->vecs_drdp) { 790c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 7910e11c21fSHong Zhang } 79213b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 793c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 7940e11c21fSHong Zhang if (ts->vecs_drdp) { 79513b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 7960e11c21fSHong Zhang } 79713b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 79813b1b0ffSHong Zhang } 79913b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 80013b1b0ffSHong Zhang } else { 801c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 8020e11c21fSHong Zhang if (ts->vecs_drdp) { 803c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 8040e11c21fSHong Zhang } 80513b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 806c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 8070e11c21fSHong Zhang if (ts->vecs_drdp) { 80813b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 8090e11c21fSHong Zhang } 81013b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 81113b1b0ffSHong Zhang } 81213b1b0ffSHong Zhang } 81313b1b0ffSHong Zhang } else { 81413b1b0ffSHong Zhang if (!th->endpoint) { 81513b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 816715f1b00SHong Zhang } 817715f1b00SHong Zhang } 8181566a47fSLisandro Dalcin PetscFunctionReturn(0); 8191566a47fSLisandro Dalcin } 8201566a47fSLisandro Dalcin 8216affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 8226affb6f8SHong Zhang { 8236affb6f8SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 8246affb6f8SHong Zhang 8256affb6f8SHong Zhang PetscFunctionBegin; 8266affb6f8SHong Zhang if (ns) *ns = 1; 8276affb6f8SHong Zhang if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 8286affb6f8SHong Zhang PetscFunctionReturn(0); 8296affb6f8SHong Zhang } 8306affb6f8SHong Zhang 831316643e7SJed Brown /*------------------------------------------------------------*/ 832277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 833316643e7SJed Brown { 834316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 835316643e7SJed Brown PetscErrorCode ierr; 836316643e7SJed Brown 837316643e7SJed Brown PetscFunctionBegin; 8386bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 8396bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 8403b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 841eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 8421566a47fSLisandro Dalcin 8431566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 8441566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 8451566a47fSLisandro Dalcin 846b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 847715f1b00SHong Zhang if (ts->forward_solve) { 848715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 849715f1b00SHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 8506f92202bSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 851715f1b00SHong Zhang } 852b5b4017aSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 85313b1b0ffSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 85413b1b0ffSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 855715f1b00SHong Zhang } 856b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 857b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 858b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 859b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 860b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 861b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 862715f1b00SHong Zhang 863277b19d0SLisandro Dalcin PetscFunctionReturn(0); 864277b19d0SLisandro Dalcin } 865277b19d0SLisandro Dalcin 866ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts) 867ece46509SHong Zhang { 868ece46509SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 869ece46509SHong Zhang PetscErrorCode ierr; 870ece46509SHong Zhang 871ece46509SHong Zhang PetscFunctionBegin; 872ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 873ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 874ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 875ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 876ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 877ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 878ece46509SHong Zhang PetscFunctionReturn(0); 879ece46509SHong Zhang } 880ece46509SHong Zhang 881277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 882277b19d0SLisandro Dalcin { 883277b19d0SLisandro Dalcin PetscErrorCode ierr; 884277b19d0SLisandro Dalcin 885277b19d0SLisandro Dalcin PetscFunctionBegin; 886277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 887b3a6b972SJed Brown if (ts->dm) { 888b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 889b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 890b3a6b972SJed Brown } 891277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 892bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 893bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 894bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 895bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 896316643e7SJed Brown PetscFunctionReturn(0); 897316643e7SJed Brown } 898316643e7SJed Brown 899316643e7SJed Brown /* 900316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9012b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 902316643e7SJed Brown */ 9030f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 904316643e7SJed Brown { 905316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 906316643e7SJed Brown PetscErrorCode ierr; 9077445fe48SJed Brown Vec X0,Xdot; 9087445fe48SJed Brown DM dm,dmsave; 9091566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 910316643e7SJed Brown 911316643e7SJed Brown PetscFunctionBegin; 9127445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 9135a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9147445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 915b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 9167445fe48SJed Brown 9177445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9187445fe48SJed Brown dmsave = ts->dm; 9197445fe48SJed Brown ts->dm = dm; 9207445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 9217445fe48SJed Brown ts->dm = dmsave; 9220d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 923316643e7SJed Brown PetscFunctionReturn(0); 924316643e7SJed Brown } 925316643e7SJed Brown 926d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 927316643e7SJed Brown { 928316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 929316643e7SJed Brown PetscErrorCode ierr; 9307445fe48SJed Brown Vec Xdot; 9317445fe48SJed Brown DM dm,dmsave; 9321566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 933316643e7SJed Brown 934316643e7SJed Brown PetscFunctionBegin; 9357445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 936be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 9370298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 9387445fe48SJed Brown 9397445fe48SJed Brown dmsave = ts->dm; 9407445fe48SJed Brown ts->dm = dm; 941d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 9427445fe48SJed Brown ts->dm = dmsave; 9430298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 944316643e7SJed Brown PetscFunctionReturn(0); 945316643e7SJed Brown } 946316643e7SJed Brown 947715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 948715f1b00SHong Zhang { 949715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 950715f1b00SHong Zhang PetscErrorCode ierr; 951715f1b00SHong Zhang 952715f1b00SHong Zhang PetscFunctionBegin; 953022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 95413b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 95513b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 956715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 957715f1b00SHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 958715f1b00SHong Zhang } 9596f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 96013b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 96113b1b0ffSHong Zhang 9626f92202bSHong Zhang if (ts->vecs_integral_sensip) { 9636f92202bSHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 9646f92202bSHong Zhang } 965b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 966b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 967715f1b00SHong Zhang PetscFunctionReturn(0); 968715f1b00SHong Zhang } 969715f1b00SHong Zhang 9707adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts) 9717adebddeSHong Zhang { 9727adebddeSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 9737adebddeSHong Zhang PetscErrorCode ierr; 9747adebddeSHong Zhang 9757adebddeSHong Zhang PetscFunctionBegin; 9767adebddeSHong Zhang if (ts->vecs_integral_sensip) { 9777adebddeSHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 9787adebddeSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 9797adebddeSHong Zhang } 9807adebddeSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 9817adebddeSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 9827adebddeSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 9837adebddeSHong Zhang PetscFunctionReturn(0); 9847adebddeSHong Zhang } 9857adebddeSHong Zhang 986316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 987316643e7SJed Brown { 988316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 9892ffb9264SLisandro Dalcin PetscBool match; 990316643e7SJed Brown PetscErrorCode ierr; 991316643e7SJed Brown 992316643e7SJed Brown PetscFunctionBegin; 993b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 994b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 995b1cb36f3SHong Zhang } 99639d13666SHong Zhang if (!th->X) { 997316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 99839d13666SHong Zhang } 99939d13666SHong Zhang if (!th->Xdot) { 1000316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 100139d13666SHong Zhang } 100239d13666SHong Zhang if (!th->X0) { 10033b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 100439d13666SHong Zhang } 10051566a47fSLisandro Dalcin if (th->endpoint) { 10061566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 10077445fe48SJed Brown } 10083b1890cdSShri Abhyankar 10091566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 10101566a47fSLisandro Dalcin 10111566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 10121566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10131566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10141566a47fSLisandro Dalcin 10151566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 10161566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 10172ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 10182ffb9264SLisandro Dalcin if (!match) { 10191566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 10201566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 10213b1890cdSShri Abhyankar } 10221566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 1023b4dd3bf9SBarry Smith PetscFunctionReturn(0); 1024b4dd3bf9SBarry Smith } 10250c7376e5SHong Zhang 1026b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1027b4dd3bf9SBarry Smith 1028b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1029b4dd3bf9SBarry Smith { 1030b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 1031b4dd3bf9SBarry Smith PetscErrorCode ierr; 1032b4dd3bf9SBarry Smith 1033b4dd3bf9SBarry Smith PetscFunctionBegin; 1034b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 1035b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 10362ca6e920SHong Zhang if (ts->vecs_sensip) { 1037b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 10382ca6e920SHong Zhang } 1039b5b4017aSHong Zhang if (ts->vecs_sensi2) { 1040b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 1041b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 1042b5b4017aSHong Zhang } 1043c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 1044c9681e5dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 1045b5b4017aSHong Zhang } 1046316643e7SJed Brown PetscFunctionReturn(0); 1047316643e7SJed Brown } 1048316643e7SJed Brown 10494416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 1050316643e7SJed Brown { 1051316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1052316643e7SJed Brown PetscErrorCode ierr; 1053316643e7SJed Brown 1054316643e7SJed Brown PetscFunctionBegin; 1055e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 1056316643e7SJed Brown { 10570298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 10580298fd71SBarry 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); 105903842d09SLisandro 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); 1060316643e7SJed Brown } 1061316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1062316643e7SJed Brown PetscFunctionReturn(0); 1063316643e7SJed Brown } 1064316643e7SJed Brown 1065316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 1066316643e7SJed Brown { 1067316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1068ace3abfcSBarry Smith PetscBool iascii; 1069316643e7SJed Brown PetscErrorCode ierr; 1070316643e7SJed Brown 1071316643e7SJed Brown PetscFunctionBegin; 1072251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1073316643e7SJed Brown if (iascii) { 10747c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 1075ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1076316643e7SJed Brown } 1077316643e7SJed Brown PetscFunctionReturn(0); 1078316643e7SJed Brown } 1079316643e7SJed Brown 1080560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 10810de4c49aSJed Brown { 10820de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 10830de4c49aSJed Brown 10840de4c49aSJed Brown PetscFunctionBegin; 10850de4c49aSJed Brown *theta = th->Theta; 10860de4c49aSJed Brown PetscFunctionReturn(0); 10870de4c49aSJed Brown } 10880de4c49aSJed Brown 1089560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 10900de4c49aSJed Brown { 10910de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 10920de4c49aSJed Brown 10930de4c49aSJed Brown PetscFunctionBegin; 10947c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 10950de4c49aSJed Brown th->Theta = theta; 10961566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 10970de4c49aSJed Brown PetscFunctionReturn(0); 10980de4c49aSJed Brown } 1099eb284becSJed Brown 1100560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 110126f2ff8fSLisandro Dalcin { 110226f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 110326f2ff8fSLisandro Dalcin 110426f2ff8fSLisandro Dalcin PetscFunctionBegin; 110526f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 110626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 110726f2ff8fSLisandro Dalcin } 110826f2ff8fSLisandro Dalcin 1109560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 1110eb284becSJed Brown { 1111eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1112eb284becSJed Brown 1113eb284becSJed Brown PetscFunctionBegin; 1114eb284becSJed Brown th->endpoint = flg; 1115eb284becSJed Brown PetscFunctionReturn(0); 1116eb284becSJed Brown } 11170de4c49aSJed Brown 1118f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1119f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 1120f9c1d6abSBarry Smith { 1121f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 1122f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 11233fd8ae06SJed Brown const PetscReal one = 1.0; 1124f9c1d6abSBarry Smith 1125f9c1d6abSBarry Smith PetscFunctionBegin; 11263fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 1127f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 1128f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 1129f9c1d6abSBarry Smith PetscFunctionReturn(0); 1130f9c1d6abSBarry Smith } 1131f9c1d6abSBarry Smith #endif 1132f9c1d6abSBarry Smith 113342682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 113442682096SHong Zhang { 113542682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 113642682096SHong Zhang 113742682096SHong Zhang PetscFunctionBegin; 11381566a47fSLisandro Dalcin if (ns) *ns = 1; 11391566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 114042682096SHong Zhang PetscFunctionReturn(0); 114142682096SHong Zhang } 1142f9c1d6abSBarry Smith 1143316643e7SJed Brown /* ------------------------------------------------------------ */ 1144316643e7SJed Brown /*MC 114596f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1146316643e7SJed Brown 1147316643e7SJed Brown Level: beginner 1148316643e7SJed Brown 11494eb428fdSBarry Smith Options Database: 1150c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1151c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 115203842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 11534eb428fdSBarry Smith 1154eb284becSJed Brown Notes: 11550c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 11560c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 11574eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 11584eb428fdSBarry Smith 1159eb284becSJed Brown This method can be applied to DAE. 1160eb284becSJed Brown 1161eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 1162eb284becSJed Brown 1163eb284becSJed Brown .vb 1164eb284becSJed Brown Theta | Theta 1165eb284becSJed Brown ------------- 1166eb284becSJed Brown | 1 1167eb284becSJed Brown .ve 1168eb284becSJed Brown 1169eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1170eb284becSJed Brown 1171eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1172eb284becSJed Brown 1173eb284becSJed Brown .vb 1174eb284becSJed Brown 0 | 0 0 1175eb284becSJed Brown 1 | 1-Theta Theta 1176eb284becSJed Brown ------------------- 1177eb284becSJed Brown | 1-Theta Theta 1178eb284becSJed Brown .ve 1179eb284becSJed Brown 1180eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1181eb284becSJed Brown 1182eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1183eb284becSJed Brown 1184eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1185eb284becSJed Brown 11864eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1187eb284becSJed Brown 1188eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1189316643e7SJed Brown 1190316643e7SJed Brown M*/ 11918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1192316643e7SJed Brown { 1193316643e7SJed Brown TS_Theta *th; 1194316643e7SJed Brown PetscErrorCode ierr; 1195316643e7SJed Brown 1196316643e7SJed Brown PetscFunctionBegin; 1197277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1198ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1199316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1200316643e7SJed Brown ts->ops->view = TSView_Theta; 1201316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 120242f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1203ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1204316643e7SJed Brown ts->ops->step = TSStep_Theta; 1205cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12061566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 120724655328SShri ts->ops->rollback = TSRollBack_Theta; 1208316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12090f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12100f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1211f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1212f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1213f9c1d6abSBarry Smith #endif 121442682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 121542f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1216b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1217b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12182ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12196affb6f8SHong Zhang 1220715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12217adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1222715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12236affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1224316643e7SJed Brown 1225efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1226efd4aadfSBarry Smith 1227b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1228316643e7SJed Brown ts->data = (void*)th; 1229316643e7SJed Brown 1230715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1231715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1232715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1233b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1234715f1b00SHong Zhang 12356f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1236316643e7SJed Brown th->Theta = 0.5; 12371566a47fSLisandro Dalcin th->order = 2; 1238bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1239bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1240bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1241bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1242316643e7SJed Brown PetscFunctionReturn(0); 1243316643e7SJed Brown } 12440de4c49aSJed Brown 12450de4c49aSJed Brown /*@ 12460de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 12470de4c49aSJed Brown 12480de4c49aSJed Brown Not Collective 12490de4c49aSJed Brown 12500de4c49aSJed Brown Input Parameter: 12510de4c49aSJed Brown . ts - timestepping context 12520de4c49aSJed Brown 12530de4c49aSJed Brown Output Parameter: 12540de4c49aSJed Brown . theta - stage abscissa 12550de4c49aSJed Brown 12560de4c49aSJed Brown Note: 12570de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 12580de4c49aSJed Brown 12590de4c49aSJed Brown Level: Advanced 12600de4c49aSJed Brown 12610de4c49aSJed Brown .seealso: TSThetaSetTheta() 12620de4c49aSJed Brown @*/ 12637087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 12640de4c49aSJed Brown { 12654ac538c5SBarry Smith PetscErrorCode ierr; 12660de4c49aSJed Brown 12670de4c49aSJed Brown PetscFunctionBegin; 1268afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12690de4c49aSJed Brown PetscValidPointer(theta,2); 12704ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 12710de4c49aSJed Brown PetscFunctionReturn(0); 12720de4c49aSJed Brown } 12730de4c49aSJed Brown 12740de4c49aSJed Brown /*@ 12750de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 12760de4c49aSJed Brown 12770de4c49aSJed Brown Not Collective 12780de4c49aSJed Brown 12790de4c49aSJed Brown Input Parameter: 12800de4c49aSJed Brown + ts - timestepping context 12810de4c49aSJed Brown - theta - stage abscissa 12820de4c49aSJed Brown 12830de4c49aSJed Brown Options Database: 12840de4c49aSJed Brown . -ts_theta_theta <theta> 12850de4c49aSJed Brown 12860de4c49aSJed Brown Level: Intermediate 12870de4c49aSJed Brown 12880de4c49aSJed Brown .seealso: TSThetaGetTheta() 12890de4c49aSJed Brown @*/ 12907087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 12910de4c49aSJed Brown { 12924ac538c5SBarry Smith PetscErrorCode ierr; 12930de4c49aSJed Brown 12940de4c49aSJed Brown PetscFunctionBegin; 1295afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12964ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 12970de4c49aSJed Brown PetscFunctionReturn(0); 12980de4c49aSJed Brown } 1299f33bbcb6SJed Brown 130026f2ff8fSLisandro Dalcin /*@ 130126f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 130226f2ff8fSLisandro Dalcin 130326f2ff8fSLisandro Dalcin Not Collective 130426f2ff8fSLisandro Dalcin 130526f2ff8fSLisandro Dalcin Input Parameter: 130626f2ff8fSLisandro Dalcin . ts - timestepping context 130726f2ff8fSLisandro Dalcin 130826f2ff8fSLisandro Dalcin Output Parameter: 130926f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 131026f2ff8fSLisandro Dalcin 131126f2ff8fSLisandro Dalcin Level: Advanced 131226f2ff8fSLisandro Dalcin 131326f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 131426f2ff8fSLisandro Dalcin @*/ 131526f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 131626f2ff8fSLisandro Dalcin { 131726f2ff8fSLisandro Dalcin PetscErrorCode ierr; 131826f2ff8fSLisandro Dalcin 131926f2ff8fSLisandro Dalcin PetscFunctionBegin; 132026f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 132126f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1322163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 132326f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 132426f2ff8fSLisandro Dalcin } 132526f2ff8fSLisandro Dalcin 1326eb284becSJed Brown /*@ 1327eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1328eb284becSJed Brown 1329eb284becSJed Brown Not Collective 1330eb284becSJed Brown 1331eb284becSJed Brown Input Parameter: 1332eb284becSJed Brown + ts - timestepping context 1333eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1334eb284becSJed Brown 1335eb284becSJed Brown Options Database: 1336eb284becSJed Brown . -ts_theta_endpoint <flg> 1337eb284becSJed Brown 1338eb284becSJed Brown Level: Intermediate 1339eb284becSJed Brown 1340eb284becSJed Brown .seealso: TSTHETA, TSCN 1341eb284becSJed Brown @*/ 1342eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1343eb284becSJed Brown { 1344eb284becSJed Brown PetscErrorCode ierr; 1345eb284becSJed Brown 1346eb284becSJed Brown PetscFunctionBegin; 1347eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1348eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1349eb284becSJed Brown PetscFunctionReturn(0); 1350eb284becSJed Brown } 1351eb284becSJed Brown 1352f33bbcb6SJed Brown /* 1353f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1354f33bbcb6SJed Brown * The creation functions for these specializations are below. 1355f33bbcb6SJed Brown */ 1356f33bbcb6SJed Brown 13571566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 13581566a47fSLisandro Dalcin { 13591566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 13601566a47fSLisandro Dalcin PetscErrorCode ierr; 13611566a47fSLisandro Dalcin 13621566a47fSLisandro Dalcin PetscFunctionBegin; 13631566a47fSLisandro 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"); 13641566a47fSLisandro 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"); 13651566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 13661566a47fSLisandro Dalcin PetscFunctionReturn(0); 13671566a47fSLisandro Dalcin } 13681566a47fSLisandro Dalcin 1369f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1370f33bbcb6SJed Brown { 1371f33bbcb6SJed Brown PetscFunctionBegin; 1372f33bbcb6SJed Brown PetscFunctionReturn(0); 1373f33bbcb6SJed Brown } 1374f33bbcb6SJed Brown 1375f33bbcb6SJed Brown /*MC 1376f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1377f33bbcb6SJed Brown 1378f33bbcb6SJed Brown Level: beginner 1379f33bbcb6SJed Brown 13804eb428fdSBarry Smith Notes: 1381c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 13824eb428fdSBarry Smith 13831566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 13844eb428fdSBarry Smith 1385f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1386f33bbcb6SJed Brown 1387f33bbcb6SJed Brown M*/ 13888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1389f33bbcb6SJed Brown { 1390f33bbcb6SJed Brown PetscErrorCode ierr; 1391f33bbcb6SJed Brown 1392f33bbcb6SJed Brown PetscFunctionBegin; 1393f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1394f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 13951566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 13960c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1397f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1398f33bbcb6SJed Brown PetscFunctionReturn(0); 1399f33bbcb6SJed Brown } 1400f33bbcb6SJed Brown 14011566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 14021566a47fSLisandro Dalcin { 14031566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14041566a47fSLisandro Dalcin PetscErrorCode ierr; 14051566a47fSLisandro Dalcin 14061566a47fSLisandro Dalcin PetscFunctionBegin; 14071566a47fSLisandro 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"); 14081566a47fSLisandro 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"); 14091566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14101566a47fSLisandro Dalcin PetscFunctionReturn(0); 14111566a47fSLisandro Dalcin } 14121566a47fSLisandro Dalcin 1413f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1414f33bbcb6SJed Brown { 1415f33bbcb6SJed Brown PetscFunctionBegin; 1416f33bbcb6SJed Brown PetscFunctionReturn(0); 1417f33bbcb6SJed Brown } 1418f33bbcb6SJed Brown 1419f33bbcb6SJed Brown /*MC 1420f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1421f33bbcb6SJed Brown 1422f33bbcb6SJed Brown Level: beginner 1423f33bbcb6SJed Brown 1424f33bbcb6SJed Brown Notes: 14257cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 14267cf5af47SJed Brown 14277cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1428f33bbcb6SJed Brown 1429f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1430f33bbcb6SJed Brown 1431f33bbcb6SJed Brown M*/ 14328cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1433f33bbcb6SJed Brown { 1434f33bbcb6SJed Brown PetscErrorCode ierr; 1435f33bbcb6SJed Brown 1436f33bbcb6SJed Brown PetscFunctionBegin; 1437f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1438f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1439eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 14400c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1441f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1442f33bbcb6SJed Brown PetscFunctionReturn(0); 1443f33bbcb6SJed Brown } 1444