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 257*c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 258*c9681e5dSHong Zhang { 259*c9681e5dSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 260*c9681e5dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 261*c9681e5dSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 262*c9681e5dSHong Zhang PetscInt nadj; 263*c9681e5dSHong Zhang Mat J,Jp; 264*c9681e5dSHong Zhang KSP ksp; 265*c9681e5dSHong Zhang PetscReal shift; 266*c9681e5dSHong Zhang PetscScalar *xarr; 267*c9681e5dSHong Zhang PetscErrorCode ierr; 268*c9681e5dSHong Zhang 269*c9681e5dSHong Zhang PetscFunctionBegin; 270*c9681e5dSHong Zhang th->status = TS_STEP_INCOMPLETE; 271*c9681e5dSHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 272*c9681e5dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 273*c9681e5dSHong Zhang 274*c9681e5dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 275*c9681e5dSHong Zhang th->stage_time = ts->ptime; /* time_step is negative*/ 276*c9681e5dSHong Zhang th->ptime = ts->ptime + ts->time_step; 277*c9681e5dSHong Zhang th->time_step = -ts->time_step; 278*c9681e5dSHong Zhang 279*c9681e5dSHong Zhang /* Build RHS for first-order adjoint */ 280*c9681e5dSHong Zhang if (ts->vec_costintegral) { /* Cost function has an integral term */ 281*c9681e5dSHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 282*c9681e5dSHong Zhang } 283*c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 284*c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 285*c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */ 286*c9681e5dSHong Zhang if (ts->vec_costintegral) { 287*c9681e5dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 288*c9681e5dSHong Zhang } 289*c9681e5dSHong Zhang } 290*c9681e5dSHong Zhang 291*c9681e5dSHong Zhang /* Build LHS for first-order adjoint */ 292*c9681e5dSHong Zhang shift = 1./th->time_step; 293*c9681e5dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 294*c9681e5dSHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 295*c9681e5dSHong Zhang 296*c9681e5dSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 297*c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 298*c9681e5dSHong Zhang KSPConvergedReason kspreason; 299*c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 300*c9681e5dSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 301*c9681e5dSHong Zhang if (kspreason < 0) { 302*c9681e5dSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 303*c9681e5dSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 304*c9681e5dSHong Zhang } 305*c9681e5dSHong Zhang } 306*c9681e5dSHong Zhang 307*c9681e5dSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 308*c9681e5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 309*c9681e5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 310*c9681e5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 311*c9681e5dSHong Zhang /* lambda_s^T F_UU w_1 */ 312*c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 313*c9681e5dSHong Zhang if (ts->vecs_fup) { 314*c9681e5dSHong Zhang /* lambda_s^T F_UP w_2 */ 315*c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 316*c9681e5dSHong Zhang } 317*c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 318*c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 319*c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 320*c9681e5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 321*c9681e5dSHong Zhang if (ts->vecs_fup) { 322*c9681e5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 323*c9681e5dSHong Zhang } 324*c9681e5dSHong Zhang } 325*c9681e5dSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 326*c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 327*c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 328*c9681e5dSHong Zhang } 329*c9681e5dSHong Zhang } 330*c9681e5dSHong Zhang 331*c9681e5dSHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 332*c9681e5dSHong Zhang shift = 0.0; 333*c9681e5dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 334*c9681e5dSHong Zhang ierr = MatScale(J,-1.);CHKERRQ(ierr); 335*c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 336*c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 337*c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr); 338*c9681e5dSHong Zhang ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 339*c9681e5dSHong Zhang if (ts->vecs_sensi2) { 340*c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 341*c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr); 342*c9681e5dSHong Zhang ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 343*c9681e5dSHong Zhang } 344*c9681e5dSHong Zhang } 345*c9681e5dSHong Zhang if (ts->vecs_sensip) { 346*c9681e5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */ 347*c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 348*c9681e5dSHong Zhang if (ts->vecs_fpu) { 349*c9681e5dSHong Zhang /* lambda_s^T F_PU w_1 */ 350*c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 351*c9681e5dSHong Zhang } 352*c9681e5dSHong Zhang if (ts->vecs_fpp) { 353*c9681e5dSHong Zhang /* lambda_s^T F_PU w_2 */ 354*c9681e5dSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 355*c9681e5dSHong Zhang } 356*c9681e5dSHong Zhang } 357*c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 358*c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 359*c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 360*c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 361*c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 362*c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 363*c9681e5dSHong Zhang if (ts->vecs_fpu) { 364*c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 365*c9681e5dSHong Zhang } 366*c9681e5dSHong Zhang if (ts->vecs_fpp) { 367*c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 368*c9681e5dSHong Zhang } 369*c9681e5dSHong Zhang } 370*c9681e5dSHong Zhang } 371*c9681e5dSHong Zhang if (ts->vec_costintegral) { 372*c9681e5dSHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 373*c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 374*c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 375*c9681e5dSHong Zhang } 376*c9681e5dSHong Zhang } 377*c9681e5dSHong Zhang } 378*c9681e5dSHong Zhang 379*c9681e5dSHong Zhang if (ts->vecs_sensi2) { 380*c9681e5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 381*c9681e5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 382*c9681e5dSHong Zhang } 383*c9681e5dSHong Zhang th->status = TS_STEP_COMPLETE; 384*c9681e5dSHong Zhang PetscFunctionReturn(0); 385*c9681e5dSHong Zhang } 386*c9681e5dSHong Zhang 38742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 3882ca6e920SHong Zhang { 3892ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 390b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 391b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 3922ca6e920SHong Zhang PetscInt nadj; 3932ca6e920SHong Zhang Mat J,Jp; 3942ca6e920SHong Zhang KSP ksp; 3952ca6e920SHong Zhang PetscReal shift; 396b5b4017aSHong Zhang PetscScalar *xarr; 397b5b4017aSHong Zhang PetscErrorCode ierr; 3982ca6e920SHong Zhang 3992ca6e920SHong Zhang PetscFunctionBegin; 400*c9681e5dSHong Zhang if (th->Theta == 1.) { 401*c9681e5dSHong Zhang ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr); 402*c9681e5dSHong Zhang PetscFunctionReturn(0); 403*c9681e5dSHong Zhang } 4042ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 405302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 4062ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 407b5e0532dSHong Zhang 408b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 409022c081aSHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 410b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 411022c081aSHong Zhang th->time_step = -ts->time_step; 4122ca6e920SHong Zhang 413b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 4142c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 41505755b9cSHong Zhang if (th->endpoint) { 416c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 41705755b9cSHong Zhang } else { 418c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 41905755b9cSHong Zhang } 42005755b9cSHong Zhang } 421abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 422b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 423022c081aSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr); 4242c39e106SBarry Smith if (ts->vec_costintegral) { 425c9ad9de0SHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 42636eaed60SHong Zhang } 4272ca6e920SHong Zhang } 4283c54f92cSHong Zhang 429b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 430022c081aSHong Zhang shift = 1./(th->Theta*th->time_step); 4313c54f92cSHong Zhang if (th->endpoint) { 432022c081aSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 4333c54f92cSHong Zhang } else { 4343c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 4353c54f92cSHong Zhang } 4362ca6e920SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 4372ca6e920SHong Zhang 438b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 439abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 440b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 4412ca6e920SHong Zhang } 4423c54f92cSHong Zhang 443b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 444b5b4017aSHong Zhang if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 445b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 446b5b4017aSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 447b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 448b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 449b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 450b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 451b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 452b5b4017aSHong Zhang if (ts->vecs_fup) { 453b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 454b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 455b5b4017aSHong Zhang } 456b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 457b5b4017aSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 458b5b4017aSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],1./shift);CHKERRQ(ierr); 459b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 460b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 461b5b4017aSHong Zhang if (ts->vecs_fup) { 462b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 463b5b4017aSHong Zhang } 464b5b4017aSHong Zhang if (ts->vec_costintegral) { 465b5b4017aSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 466b5b4017aSHong Zhang } 467b5b4017aSHong Zhang } 468b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 469b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 470b5b4017aSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 471b5b4017aSHong Zhang } 472b5b4017aSHong Zhang } 473b5b4017aSHong Zhang 47436eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 4756affb6f8SHong Zhang if(th->endpoint) { /* two-stage Theta methods */ 476b5b4017aSHong Zhang if (th->Theta!=1.) { /* general case */ 477022c081aSHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 478b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 479b5b4017aSHong Zhang if (ts->vec_costintegral) { /* R_U at t_n */ 480c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 4818263dbbfSHong Zhang } 482abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 483b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 4843c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 485b5b4017aSHong Zhang if (ts->vec_costintegral) { 486b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1./shift,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 487b5b4017aSHong Zhang } 488b5b4017aSHong Zhang } 489b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* second-order */ 490b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 491b5b4017aSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 492b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 493b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 494b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 495b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 496b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 497b5b4017aSHong Zhang if (ts->vecs_fup) { 498b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 499b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 500b5b4017aSHong Zhang } 501b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 502b5b4017aSHong Zhang /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) R_U */ 503b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 504b5b4017aSHong Zhang ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr); 505b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 506b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 507b5b4017aSHong Zhang if (ts->vecs_fup) { 508b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fup[nadj]);CHKERRQ(ierr); 509b5b4017aSHong Zhang } 510b5b4017aSHong Zhang if (ts->vec_costintegral) { 511b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 512b5b4017aSHong Zhang } 513b5b4017aSHong Zhang } 5142ca6e920SHong Zhang } 515b5e0532dSHong Zhang } else { /* backward Euler */ 516b5e0532dSHong Zhang shift = 0.0; 517b5b4017aSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_u */ 518b5e0532dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 519b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 520022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 521b5b4017aSHong Zhang if (ts->vec_costintegral) { /* wrong? */ 522c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 523b5e0532dSHong Zhang } 524b5e0532dSHong Zhang } 525b5b4017aSHong Zhang if (ts->vecs_sensi2) { 526b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 527b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 528b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],-th->time_step,VecsSensi2Temp[nadj]);CHKERRQ(ierr); 529b5b4017aSHong Zhang } 530b5b4017aSHong Zhang } 531b5e0532dSHong Zhang } 5323fd52205SHong Zhang 5333fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 534b5b4017aSHong Zhang /* U_{n+1} */ 535c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 536b5b4017aSHong Zhang if (ts->vec_costintegral) { 537b5b4017aSHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 538b5b4017aSHong Zhang } 539abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 540b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 541022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 5423fd52205SHong Zhang } 543b5b4017aSHong Zhang if (ts->vecs_sensip2) { /* second-order */ 544b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 545b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 546b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 547b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 548b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 549b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 550b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 551b5b4017aSHong Zhang if (ts->vecs_fpu) { 552b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 553b5b4017aSHong Zhang } 554b5b4017aSHong Zhang if (ts->vecs_fpp) { 555b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 556b5b4017aSHong Zhang } 557b5b4017aSHong Zhang if (ts->vec_costintegral) { 558b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 559b5b4017aSHong Zhang } 560b5b4017aSHong Zhang } 561b5b4017aSHong Zhang } 562b5b4017aSHong Zhang 563b5b4017aSHong Zhang /* U_s */ 564b5e0532dSHong Zhang if (th->Theta!=1.) { 565c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 566b5b4017aSHong Zhang if (ts->vec_costintegral) { 567b5b4017aSHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 568b5b4017aSHong Zhang } 569abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 570b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 571022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 572b5b4017aSHong Zhang if (ts->vecs_sensip2) { /* second-order */ 573b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 574b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 575b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 576b5b4017aSHong Zhang ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 577b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 578b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 579b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 580b5b4017aSHong Zhang if (ts->vecs_fpu) { 581b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 582b5e0532dSHong Zhang } 583b5b4017aSHong Zhang if (ts->vecs_fpp) { 584b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 5853fd52205SHong Zhang } 5862c39e106SBarry Smith if (ts->vec_costintegral) { 587b5b4017aSHong Zhang ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 58836eaed60SHong Zhang } 589b5b4017aSHong Zhang } 59036eaed60SHong Zhang } 59136eaed60SHong Zhang } 5923fd52205SHong Zhang } 593b5e0532dSHong Zhang } 5943fd52205SHong Zhang } else { /* one-stage case */ 5953c54f92cSHong Zhang shift = 0.0; 596a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 5972c39e106SBarry Smith if (ts->vec_costintegral) { 598c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 5993c54f92cSHong Zhang } 600abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 601b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 602022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 6032c39e106SBarry Smith if (ts->vec_costintegral) { 604c9ad9de0SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 60536eaed60SHong Zhang } 6062ca6e920SHong Zhang } 6073fd52205SHong Zhang if (ts->vecs_sensip) { 608c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 609abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 610b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 611022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 6123fd52205SHong Zhang } 6132c39e106SBarry Smith if (ts->vec_costintegral) { 614c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 615abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 616022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 61736eaed60SHong Zhang } 61836eaed60SHong Zhang } 6193fd52205SHong Zhang } 6202ca6e920SHong Zhang } 6212ca6e920SHong Zhang 6222ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6232ca6e920SHong Zhang PetscFunctionReturn(0); 6242ca6e920SHong Zhang } 6252ca6e920SHong Zhang 626cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 627cd652676SJed Brown { 628cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 629be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 630cd652676SJed Brown PetscErrorCode ierr; 631cd652676SJed Brown 632cd652676SJed Brown PetscFunctionBegin; 633a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 634be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 635be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 636cd652676SJed Brown PetscFunctionReturn(0); 637cd652676SJed Brown } 638cd652676SJed Brown 6391566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 6401566a47fSLisandro Dalcin { 6411566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 6421566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 6431566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 6447453f775SEmil Constantinescu PetscReal wltea,wlter; 6451566a47fSLisandro Dalcin PetscErrorCode ierr; 6461566a47fSLisandro Dalcin 6471566a47fSLisandro Dalcin PetscFunctionBegin; 6482ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 6491566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 650fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 6511566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 652fecfb714SLisandro Dalcin { 653be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 654be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 6551566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 6561566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 6571566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 6581566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 6591566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 6607453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 6611566a47fSLisandro Dalcin } 6621566a47fSLisandro Dalcin if (order) *order = 2; 6631566a47fSLisandro Dalcin PetscFunctionReturn(0); 6641566a47fSLisandro Dalcin } 6651566a47fSLisandro Dalcin 6661566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 6671566a47fSLisandro Dalcin { 6681566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 66913b1b0ffSHong Zhang PetscInt ncost; 6701566a47fSLisandro Dalcin PetscErrorCode ierr; 6711566a47fSLisandro Dalcin 6721566a47fSLisandro Dalcin PetscFunctionBegin; 6731566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 6741566a47fSLisandro Dalcin if (ts->vec_costintegral && ts->costintegralfwd) { 6751566a47fSLisandro Dalcin ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 6761566a47fSLisandro Dalcin } 677715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 67813b1b0ffSHong Zhang if (ts->mat_sensip) { 67913b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 6806f92202bSHong Zhang } 6816f92202bSHong Zhang if (ts->vecs_integral_sensip) { 6826f92202bSHong Zhang for (ncost=0;ncost<ts->numcost;ncost++) { 6836f92202bSHong Zhang ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 6846f92202bSHong Zhang } 6856f92202bSHong Zhang } 686715f1b00SHong Zhang PetscFunctionReturn(0); 687715f1b00SHong Zhang } 688715f1b00SHong Zhang 689715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 690715f1b00SHong Zhang { 691715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 69213b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 693b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 69413b1b0ffSHong Zhang PetscInt ncost,ntlm; 695715f1b00SHong Zhang KSP ksp; 696715f1b00SHong Zhang Mat J,Jp; 697715f1b00SHong Zhang PetscReal shift; 69813b1b0ffSHong Zhang PetscScalar *barr,*xarr; 699715f1b00SHong Zhang PetscErrorCode ierr; 700715f1b00SHong Zhang 701715f1b00SHong Zhang PetscFunctionBegin; 70213b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 70313b1b0ffSHong Zhang 7046f92202bSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 7056f92202bSHong Zhang if (ts->vecs_integral_sensip) { 7066f92202bSHong Zhang ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 7076f92202bSHong Zhang } 7086f92202bSHong Zhang } 7096f92202bSHong Zhang 710715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 711715f1b00SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 712715f1b00SHong Zhang 713715f1b00SHong Zhang /* Build RHS */ 714715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 715715f1b00SHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 716715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,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,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 719715f1b00SHong Zhang 720022c081aSHong Zhang /* Add the f_p forcing terms */ 7210e11c21fSHong Zhang if (ts->Jacp) { 722c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 72313b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(1.-th->Theta)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 724c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 72513b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 7260e11c21fSHong Zhang } 727715f1b00SHong Zhang } else { /* 1-stage method */ 728715f1b00SHong Zhang shift = 0.0; 729715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 73013b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 73113b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 73213b1b0ffSHong Zhang 733022c081aSHong Zhang /* Add the f_p forcing terms */ 7340e11c21fSHong Zhang if (ts->Jacp) { 735c9ad9de0SHong Zhang ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 73613b1b0ffSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 737715f1b00SHong Zhang } 7380e11c21fSHong Zhang } 739715f1b00SHong Zhang 740715f1b00SHong Zhang /* Build LHS */ 741715f1b00SHong Zhang shift = 1/(th->Theta*th->time_step); 742715f1b00SHong Zhang if (th->endpoint) { 743715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 744715f1b00SHong Zhang } else { 745715f1b00SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 746715f1b00SHong Zhang } 747715f1b00SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 748715f1b00SHong Zhang 749715f1b00SHong Zhang /* 750715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 751c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 752715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 753715f1b00SHong Zhang */ 754715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 755715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 756c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 7570e11c21fSHong Zhang if (ts->vecs_drdp) { 758c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 7590e11c21fSHong Zhang } 760715f1b00SHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 761c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 7620e11c21fSHong Zhang if (ts->vecs_drdp) { 76313b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 7640e11c21fSHong Zhang } 765715f1b00SHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 766715f1b00SHong Zhang } 767715f1b00SHong Zhang } 768715f1b00SHong Zhang } 769715f1b00SHong Zhang 770715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 771022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 772b5b4017aSHong Zhang KSPConvergedReason kspreason; 77313b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 774b5b4017aSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 775715f1b00SHong Zhang if (th->endpoint) { 77613b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 777b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 778b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 779b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 78013b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 781715f1b00SHong Zhang } else { 782b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 783715f1b00SHong Zhang } 784b5b4017aSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 785b5b4017aSHong Zhang if (kspreason < 0) { 786b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 787b5b4017aSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 788b5b4017aSHong Zhang } 789b5b4017aSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 79013b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 791715f1b00SHong Zhang } 792715f1b00SHong Zhang 793b5b4017aSHong Zhang 79413b1b0ffSHong Zhang /* 79513b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 796c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 79713b1b0ffSHong Zhang */ 79813b1b0ffSHong Zhang if (ts->vecs_integral_sensip) { 79913b1b0ffSHong Zhang if (!th->endpoint) { 80013b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 801c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 8020e11c21fSHong Zhang if (ts->vecs_drdp) { 803c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,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->VecIntegralSensipTemp);CHKERRQ(ierr); 81113b1b0ffSHong Zhang } 81213b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 81313b1b0ffSHong Zhang } else { 814c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 8150e11c21fSHong Zhang if (ts->vecs_drdp) { 816c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 8170e11c21fSHong Zhang } 81813b1b0ffSHong Zhang for (ncost=0; ncost<ts->numcost; ncost++) { 819c9ad9de0SHong Zhang ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 8200e11c21fSHong Zhang if (ts->vecs_drdp) { 82113b1b0ffSHong Zhang ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 8220e11c21fSHong Zhang } 82313b1b0ffSHong Zhang ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 82413b1b0ffSHong Zhang } 82513b1b0ffSHong Zhang } 82613b1b0ffSHong Zhang } else { 82713b1b0ffSHong Zhang if (!th->endpoint) { 82813b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 829715f1b00SHong Zhang } 830715f1b00SHong Zhang } 8311566a47fSLisandro Dalcin PetscFunctionReturn(0); 8321566a47fSLisandro Dalcin } 8331566a47fSLisandro Dalcin 8346affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 8356affb6f8SHong Zhang { 8366affb6f8SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 8376affb6f8SHong Zhang 8386affb6f8SHong Zhang PetscFunctionBegin; 8396affb6f8SHong Zhang if (ns) *ns = 1; 8406affb6f8SHong Zhang if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 8416affb6f8SHong Zhang PetscFunctionReturn(0); 8426affb6f8SHong Zhang } 8436affb6f8SHong Zhang 844316643e7SJed Brown /*------------------------------------------------------------*/ 845277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 846316643e7SJed Brown { 847316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 848316643e7SJed Brown PetscErrorCode ierr; 849316643e7SJed Brown 850316643e7SJed Brown PetscFunctionBegin; 8516bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 8526bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 8533b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 854eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 8551566a47fSLisandro Dalcin 8561566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 8571566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 8581566a47fSLisandro Dalcin 859b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 860715f1b00SHong Zhang if (ts->forward_solve) { 861715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 862715f1b00SHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 8636f92202bSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 864715f1b00SHong Zhang } 865b5b4017aSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 86613b1b0ffSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 86713b1b0ffSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 868715f1b00SHong Zhang } 869b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 870b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 871b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 872b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 873b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 874b5b4017aSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 875715f1b00SHong Zhang 876277b19d0SLisandro Dalcin PetscFunctionReturn(0); 877277b19d0SLisandro Dalcin } 878277b19d0SLisandro Dalcin 879ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts) 880ece46509SHong Zhang { 881ece46509SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 882ece46509SHong Zhang PetscErrorCode ierr; 883ece46509SHong Zhang 884ece46509SHong Zhang PetscFunctionBegin; 885ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 886ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 887ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 888ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 889ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 890ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 891ece46509SHong Zhang PetscFunctionReturn(0); 892ece46509SHong Zhang } 893ece46509SHong Zhang 894277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 895277b19d0SLisandro Dalcin { 896277b19d0SLisandro Dalcin PetscErrorCode ierr; 897277b19d0SLisandro Dalcin 898277b19d0SLisandro Dalcin PetscFunctionBegin; 899277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 900b3a6b972SJed Brown if (ts->dm) { 901b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 902b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 903b3a6b972SJed Brown } 904277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 905bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 906bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 907bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 908bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 909316643e7SJed Brown PetscFunctionReturn(0); 910316643e7SJed Brown } 911316643e7SJed Brown 912316643e7SJed Brown /* 913316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9142b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 915316643e7SJed Brown */ 9160f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 917316643e7SJed Brown { 918316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 919316643e7SJed Brown PetscErrorCode ierr; 9207445fe48SJed Brown Vec X0,Xdot; 9217445fe48SJed Brown DM dm,dmsave; 9221566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 923316643e7SJed Brown 924316643e7SJed Brown PetscFunctionBegin; 9257445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 9265a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9277445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 928b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 9297445fe48SJed Brown 9307445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9317445fe48SJed Brown dmsave = ts->dm; 9327445fe48SJed Brown ts->dm = dm; 9337445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 9347445fe48SJed Brown ts->dm = dmsave; 9350d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 936316643e7SJed Brown PetscFunctionReturn(0); 937316643e7SJed Brown } 938316643e7SJed Brown 939d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 940316643e7SJed Brown { 941316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 942316643e7SJed Brown PetscErrorCode ierr; 9437445fe48SJed Brown Vec Xdot; 9447445fe48SJed Brown DM dm,dmsave; 9451566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 946316643e7SJed Brown 947316643e7SJed Brown PetscFunctionBegin; 9487445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 949be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 9500298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 9517445fe48SJed Brown 9527445fe48SJed Brown dmsave = ts->dm; 9537445fe48SJed Brown ts->dm = dm; 954d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 9557445fe48SJed Brown ts->dm = dmsave; 9560298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 957316643e7SJed Brown PetscFunctionReturn(0); 958316643e7SJed Brown } 959316643e7SJed Brown 960715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 961715f1b00SHong Zhang { 962715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 963715f1b00SHong Zhang PetscErrorCode ierr; 964715f1b00SHong Zhang 965715f1b00SHong Zhang PetscFunctionBegin; 966022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 96713b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 96813b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 969715f1b00SHong Zhang if (ts->vecs_integral_sensip) { 970715f1b00SHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 971715f1b00SHong Zhang } 9726f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 97313b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 97413b1b0ffSHong Zhang 9756f92202bSHong Zhang if (ts->vecs_integral_sensip) { 9766f92202bSHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 9776f92202bSHong Zhang } 978b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 979b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 980715f1b00SHong Zhang PetscFunctionReturn(0); 981715f1b00SHong Zhang } 982715f1b00SHong Zhang 9837adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts) 9847adebddeSHong Zhang { 9857adebddeSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 9867adebddeSHong Zhang PetscErrorCode ierr; 9877adebddeSHong Zhang 9887adebddeSHong Zhang PetscFunctionBegin; 9897adebddeSHong Zhang if (ts->vecs_integral_sensip) { 9907adebddeSHong Zhang ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 9917adebddeSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 9927adebddeSHong Zhang } 9937adebddeSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 9947adebddeSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 9957adebddeSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 9967adebddeSHong Zhang PetscFunctionReturn(0); 9977adebddeSHong Zhang } 9987adebddeSHong Zhang 999316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 1000316643e7SJed Brown { 1001316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 10022ffb9264SLisandro Dalcin PetscBool match; 1003316643e7SJed Brown PetscErrorCode ierr; 1004316643e7SJed Brown 1005316643e7SJed Brown PetscFunctionBegin; 1006b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 1007b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 1008b1cb36f3SHong Zhang } 100939d13666SHong Zhang if (!th->X) { 1010316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 101139d13666SHong Zhang } 101239d13666SHong Zhang if (!th->Xdot) { 1013316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 101439d13666SHong Zhang } 101539d13666SHong Zhang if (!th->X0) { 10163b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 101739d13666SHong Zhang } 10181566a47fSLisandro Dalcin if (th->endpoint) { 10191566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 10207445fe48SJed Brown } 10213b1890cdSShri Abhyankar 10221566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 10231566a47fSLisandro Dalcin 10241566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 10251566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10261566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10271566a47fSLisandro Dalcin 10281566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 10291566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 10302ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 10312ffb9264SLisandro Dalcin if (!match) { 10321566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 10331566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 10343b1890cdSShri Abhyankar } 10351566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 1036b4dd3bf9SBarry Smith PetscFunctionReturn(0); 1037b4dd3bf9SBarry Smith } 10380c7376e5SHong Zhang 1039b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1040b4dd3bf9SBarry Smith 1041b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1042b4dd3bf9SBarry Smith { 1043b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 1044b4dd3bf9SBarry Smith PetscErrorCode ierr; 1045b4dd3bf9SBarry Smith 1046b4dd3bf9SBarry Smith PetscFunctionBegin; 1047b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 1048b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 10492ca6e920SHong Zhang if (ts->vecs_sensip) { 1050b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 10512ca6e920SHong Zhang } 1052b5b4017aSHong Zhang if (ts->vecs_sensi2) { 1053b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 1054b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 1055b5b4017aSHong Zhang } 1056*c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 1057*c9681e5dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 1058b5b4017aSHong Zhang } 1059316643e7SJed Brown PetscFunctionReturn(0); 1060316643e7SJed Brown } 1061316643e7SJed Brown 10624416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 1063316643e7SJed Brown { 1064316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1065316643e7SJed Brown PetscErrorCode ierr; 1066316643e7SJed Brown 1067316643e7SJed Brown PetscFunctionBegin; 1068e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 1069316643e7SJed Brown { 10700298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 10710298fd71SBarry 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); 107203842d09SLisandro 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); 1073316643e7SJed Brown } 1074316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1075316643e7SJed Brown PetscFunctionReturn(0); 1076316643e7SJed Brown } 1077316643e7SJed Brown 1078316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 1079316643e7SJed Brown { 1080316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1081ace3abfcSBarry Smith PetscBool iascii; 1082316643e7SJed Brown PetscErrorCode ierr; 1083316643e7SJed Brown 1084316643e7SJed Brown PetscFunctionBegin; 1085251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1086316643e7SJed Brown if (iascii) { 10877c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 1088ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1089316643e7SJed Brown } 1090316643e7SJed Brown PetscFunctionReturn(0); 1091316643e7SJed Brown } 1092316643e7SJed Brown 1093560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 10940de4c49aSJed Brown { 10950de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 10960de4c49aSJed Brown 10970de4c49aSJed Brown PetscFunctionBegin; 10980de4c49aSJed Brown *theta = th->Theta; 10990de4c49aSJed Brown PetscFunctionReturn(0); 11000de4c49aSJed Brown } 11010de4c49aSJed Brown 1102560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 11030de4c49aSJed Brown { 11040de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11050de4c49aSJed Brown 11060de4c49aSJed Brown PetscFunctionBegin; 11077c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 11080de4c49aSJed Brown th->Theta = theta; 11091566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 11100de4c49aSJed Brown PetscFunctionReturn(0); 11110de4c49aSJed Brown } 1112eb284becSJed Brown 1113560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 111426f2ff8fSLisandro Dalcin { 111526f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 111626f2ff8fSLisandro Dalcin 111726f2ff8fSLisandro Dalcin PetscFunctionBegin; 111826f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 111926f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 112026f2ff8fSLisandro Dalcin } 112126f2ff8fSLisandro Dalcin 1122560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 1123eb284becSJed Brown { 1124eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1125eb284becSJed Brown 1126eb284becSJed Brown PetscFunctionBegin; 1127eb284becSJed Brown th->endpoint = flg; 1128eb284becSJed Brown PetscFunctionReturn(0); 1129eb284becSJed Brown } 11300de4c49aSJed Brown 1131f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1132f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 1133f9c1d6abSBarry Smith { 1134f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 1135f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 11363fd8ae06SJed Brown const PetscReal one = 1.0; 1137f9c1d6abSBarry Smith 1138f9c1d6abSBarry Smith PetscFunctionBegin; 11393fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 1140f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 1141f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 1142f9c1d6abSBarry Smith PetscFunctionReturn(0); 1143f9c1d6abSBarry Smith } 1144f9c1d6abSBarry Smith #endif 1145f9c1d6abSBarry Smith 114642682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 114742682096SHong Zhang { 114842682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 114942682096SHong Zhang 115042682096SHong Zhang PetscFunctionBegin; 11511566a47fSLisandro Dalcin if (ns) *ns = 1; 11521566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 115342682096SHong Zhang PetscFunctionReturn(0); 115442682096SHong Zhang } 1155f9c1d6abSBarry Smith 1156316643e7SJed Brown /* ------------------------------------------------------------ */ 1157316643e7SJed Brown /*MC 115896f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1159316643e7SJed Brown 1160316643e7SJed Brown Level: beginner 1161316643e7SJed Brown 11624eb428fdSBarry Smith Options Database: 1163c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1164c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 116503842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 11664eb428fdSBarry Smith 1167eb284becSJed Brown Notes: 11680c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 11690c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 11704eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 11714eb428fdSBarry Smith 1172eb284becSJed Brown This method can be applied to DAE. 1173eb284becSJed Brown 1174eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 1175eb284becSJed Brown 1176eb284becSJed Brown .vb 1177eb284becSJed Brown Theta | Theta 1178eb284becSJed Brown ------------- 1179eb284becSJed Brown | 1 1180eb284becSJed Brown .ve 1181eb284becSJed Brown 1182eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1183eb284becSJed Brown 1184eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1185eb284becSJed Brown 1186eb284becSJed Brown .vb 1187eb284becSJed Brown 0 | 0 0 1188eb284becSJed Brown 1 | 1-Theta Theta 1189eb284becSJed Brown ------------------- 1190eb284becSJed Brown | 1-Theta Theta 1191eb284becSJed Brown .ve 1192eb284becSJed Brown 1193eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1194eb284becSJed Brown 1195eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1196eb284becSJed Brown 1197eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1198eb284becSJed Brown 11994eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1200eb284becSJed Brown 1201eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1202316643e7SJed Brown 1203316643e7SJed Brown M*/ 12048cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1205316643e7SJed Brown { 1206316643e7SJed Brown TS_Theta *th; 1207316643e7SJed Brown PetscErrorCode ierr; 1208316643e7SJed Brown 1209316643e7SJed Brown PetscFunctionBegin; 1210277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1211ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1212316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1213316643e7SJed Brown ts->ops->view = TSView_Theta; 1214316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 121542f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1216ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1217316643e7SJed Brown ts->ops->step = TSStep_Theta; 1218cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12191566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 122024655328SShri ts->ops->rollback = TSRollBack_Theta; 1221316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12220f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12230f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1224f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1225f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1226f9c1d6abSBarry Smith #endif 122742682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 122842f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1229b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1230b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12312ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12326affb6f8SHong Zhang 1233715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12347adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1235715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12366affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1237316643e7SJed Brown 1238efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1239efd4aadfSBarry Smith 1240b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1241316643e7SJed Brown ts->data = (void*)th; 1242316643e7SJed Brown 1243715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1244715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1245715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1246b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1247715f1b00SHong Zhang 12486f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1249316643e7SJed Brown th->Theta = 0.5; 12501566a47fSLisandro Dalcin th->order = 2; 1251bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1252bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1253bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1254bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1255316643e7SJed Brown PetscFunctionReturn(0); 1256316643e7SJed Brown } 12570de4c49aSJed Brown 12580de4c49aSJed Brown /*@ 12590de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 12600de4c49aSJed Brown 12610de4c49aSJed Brown Not Collective 12620de4c49aSJed Brown 12630de4c49aSJed Brown Input Parameter: 12640de4c49aSJed Brown . ts - timestepping context 12650de4c49aSJed Brown 12660de4c49aSJed Brown Output Parameter: 12670de4c49aSJed Brown . theta - stage abscissa 12680de4c49aSJed Brown 12690de4c49aSJed Brown Note: 12700de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 12710de4c49aSJed Brown 12720de4c49aSJed Brown Level: Advanced 12730de4c49aSJed Brown 12740de4c49aSJed Brown .seealso: TSThetaSetTheta() 12750de4c49aSJed Brown @*/ 12767087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 12770de4c49aSJed Brown { 12784ac538c5SBarry Smith PetscErrorCode ierr; 12790de4c49aSJed Brown 12800de4c49aSJed Brown PetscFunctionBegin; 1281afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12820de4c49aSJed Brown PetscValidPointer(theta,2); 12834ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 12840de4c49aSJed Brown PetscFunctionReturn(0); 12850de4c49aSJed Brown } 12860de4c49aSJed Brown 12870de4c49aSJed Brown /*@ 12880de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 12890de4c49aSJed Brown 12900de4c49aSJed Brown Not Collective 12910de4c49aSJed Brown 12920de4c49aSJed Brown Input Parameter: 12930de4c49aSJed Brown + ts - timestepping context 12940de4c49aSJed Brown - theta - stage abscissa 12950de4c49aSJed Brown 12960de4c49aSJed Brown Options Database: 12970de4c49aSJed Brown . -ts_theta_theta <theta> 12980de4c49aSJed Brown 12990de4c49aSJed Brown Level: Intermediate 13000de4c49aSJed Brown 13010de4c49aSJed Brown .seealso: TSThetaGetTheta() 13020de4c49aSJed Brown @*/ 13037087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 13040de4c49aSJed Brown { 13054ac538c5SBarry Smith PetscErrorCode ierr; 13060de4c49aSJed Brown 13070de4c49aSJed Brown PetscFunctionBegin; 1308afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13094ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 13100de4c49aSJed Brown PetscFunctionReturn(0); 13110de4c49aSJed Brown } 1312f33bbcb6SJed Brown 131326f2ff8fSLisandro Dalcin /*@ 131426f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 131526f2ff8fSLisandro Dalcin 131626f2ff8fSLisandro Dalcin Not Collective 131726f2ff8fSLisandro Dalcin 131826f2ff8fSLisandro Dalcin Input Parameter: 131926f2ff8fSLisandro Dalcin . ts - timestepping context 132026f2ff8fSLisandro Dalcin 132126f2ff8fSLisandro Dalcin Output Parameter: 132226f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 132326f2ff8fSLisandro Dalcin 132426f2ff8fSLisandro Dalcin Level: Advanced 132526f2ff8fSLisandro Dalcin 132626f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 132726f2ff8fSLisandro Dalcin @*/ 132826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 132926f2ff8fSLisandro Dalcin { 133026f2ff8fSLisandro Dalcin PetscErrorCode ierr; 133126f2ff8fSLisandro Dalcin 133226f2ff8fSLisandro Dalcin PetscFunctionBegin; 133326f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 133426f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1335163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 133626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 133726f2ff8fSLisandro Dalcin } 133826f2ff8fSLisandro Dalcin 1339eb284becSJed Brown /*@ 1340eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1341eb284becSJed Brown 1342eb284becSJed Brown Not Collective 1343eb284becSJed Brown 1344eb284becSJed Brown Input Parameter: 1345eb284becSJed Brown + ts - timestepping context 1346eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1347eb284becSJed Brown 1348eb284becSJed Brown Options Database: 1349eb284becSJed Brown . -ts_theta_endpoint <flg> 1350eb284becSJed Brown 1351eb284becSJed Brown Level: Intermediate 1352eb284becSJed Brown 1353eb284becSJed Brown .seealso: TSTHETA, TSCN 1354eb284becSJed Brown @*/ 1355eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1356eb284becSJed Brown { 1357eb284becSJed Brown PetscErrorCode ierr; 1358eb284becSJed Brown 1359eb284becSJed Brown PetscFunctionBegin; 1360eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1361eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1362eb284becSJed Brown PetscFunctionReturn(0); 1363eb284becSJed Brown } 1364eb284becSJed Brown 1365f33bbcb6SJed Brown /* 1366f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1367f33bbcb6SJed Brown * The creation functions for these specializations are below. 1368f33bbcb6SJed Brown */ 1369f33bbcb6SJed Brown 13701566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 13711566a47fSLisandro Dalcin { 13721566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 13731566a47fSLisandro Dalcin PetscErrorCode ierr; 13741566a47fSLisandro Dalcin 13751566a47fSLisandro Dalcin PetscFunctionBegin; 13761566a47fSLisandro 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"); 13771566a47fSLisandro 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"); 13781566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 13791566a47fSLisandro Dalcin PetscFunctionReturn(0); 13801566a47fSLisandro Dalcin } 13811566a47fSLisandro Dalcin 1382f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1383f33bbcb6SJed Brown { 1384f33bbcb6SJed Brown PetscFunctionBegin; 1385f33bbcb6SJed Brown PetscFunctionReturn(0); 1386f33bbcb6SJed Brown } 1387f33bbcb6SJed Brown 1388f33bbcb6SJed Brown /*MC 1389f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1390f33bbcb6SJed Brown 1391f33bbcb6SJed Brown Level: beginner 1392f33bbcb6SJed Brown 13934eb428fdSBarry Smith Notes: 1394c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 13954eb428fdSBarry Smith 13961566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 13974eb428fdSBarry Smith 1398f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1399f33bbcb6SJed Brown 1400f33bbcb6SJed Brown M*/ 14018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1402f33bbcb6SJed Brown { 1403f33bbcb6SJed Brown PetscErrorCode ierr; 1404f33bbcb6SJed Brown 1405f33bbcb6SJed Brown PetscFunctionBegin; 1406f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1407f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 14081566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 14090c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1410f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1411f33bbcb6SJed Brown PetscFunctionReturn(0); 1412f33bbcb6SJed Brown } 1413f33bbcb6SJed Brown 14141566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 14151566a47fSLisandro Dalcin { 14161566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14171566a47fSLisandro Dalcin PetscErrorCode ierr; 14181566a47fSLisandro Dalcin 14191566a47fSLisandro Dalcin PetscFunctionBegin; 14201566a47fSLisandro 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"); 14211566a47fSLisandro 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"); 14221566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14231566a47fSLisandro Dalcin PetscFunctionReturn(0); 14241566a47fSLisandro Dalcin } 14251566a47fSLisandro Dalcin 1426f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1427f33bbcb6SJed Brown { 1428f33bbcb6SJed Brown PetscFunctionBegin; 1429f33bbcb6SJed Brown PetscFunctionReturn(0); 1430f33bbcb6SJed Brown } 1431f33bbcb6SJed Brown 1432f33bbcb6SJed Brown /*MC 1433f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1434f33bbcb6SJed Brown 1435f33bbcb6SJed Brown Level: beginner 1436f33bbcb6SJed Brown 1437f33bbcb6SJed Brown Notes: 14387cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 14397cf5af47SJed Brown 14407cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1441f33bbcb6SJed Brown 1442f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1443f33bbcb6SJed Brown 1444f33bbcb6SJed Brown M*/ 14458cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1446f33bbcb6SJed Brown { 1447f33bbcb6SJed Brown PetscErrorCode ierr; 1448f33bbcb6SJed Brown 1449f33bbcb6SJed Brown PetscFunctionBegin; 1450f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1451f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1452eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 14530c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1454f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1455f33bbcb6SJed Brown PetscFunctionReturn(0); 1456f33bbcb6SJed Brown } 1457