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; 17*3e05ceb1SHong Zhang PetscReal shift; 181566a47fSLisandro Dalcin PetscInt order; 191566a47fSLisandro Dalcin PetscBool endpoint; 201566a47fSLisandro Dalcin PetscBool extrapolate; 21715f1b00SHong Zhang TSStepStatus status; 22715f1b00SHong Zhang Vec VecCostIntegral0; /* Backup for roll-backs due to events */ 231566a47fSLisandro Dalcin 24715f1b00SHong Zhang /* context for sensitivity analysis */ 25022c081aSHong Zhang PetscInt num_tlm; /* Total number of tangent linear equations */ 26b5e0532dSHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 27b5e0532dSHong Zhang Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 2813b1b0ffSHong Zhang Vec *VecsSensiTemp; /* Vector to be multiplied with Jacobian transpose */ 2913b1b0ffSHong Zhang Mat MatDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 30b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */ 3113b1b0ffSHong Zhang Mat MatFwdSensip0; /* backup for roll-backs due to events */ 327207e0fdSHong Zhang Mat MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 337207e0fdSHong Zhang Mat MatIntegralSensip0; /* backup for roll-backs due to events */ 34b5b4017aSHong Zhang Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */ 35b5b4017aSHong Zhang Vec *VecsDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */ 36b5b4017aSHong Zhang Vec *VecsSensi2Temp; /* Working vectors that holds the residual for the second-order adjoint */ 37b5b4017aSHong Zhang Vec *VecsAffine; /* Working vectors to store residuals */ 38715f1b00SHong Zhang /* context for error estimation */ 391566a47fSLisandro Dalcin Vec vec_sol_prev; 401566a47fSLisandro Dalcin Vec vec_lte_work; 41316643e7SJed Brown } TS_Theta; 42316643e7SJed Brown 437445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 447445fe48SJed Brown { 457445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 467445fe48SJed Brown PetscErrorCode ierr; 477445fe48SJed Brown 487445fe48SJed Brown PetscFunctionBegin; 497445fe48SJed Brown if (X0) { 507445fe48SJed Brown if (dm && dm != ts->dm) { 510d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 527445fe48SJed Brown } else *X0 = ts->vec_sol; 537445fe48SJed Brown } 547445fe48SJed Brown if (Xdot) { 557445fe48SJed Brown if (dm && dm != ts->dm) { 560d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 577445fe48SJed Brown } else *Xdot = th->Xdot; 587445fe48SJed Brown } 597445fe48SJed Brown PetscFunctionReturn(0); 607445fe48SJed Brown } 617445fe48SJed Brown 620d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 630d0b770aSPeter Brune { 640d0b770aSPeter Brune PetscErrorCode ierr; 650d0b770aSPeter Brune 660d0b770aSPeter Brune PetscFunctionBegin; 670d0b770aSPeter Brune if (X0) { 680d0b770aSPeter Brune if (dm && dm != ts->dm) { 690d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 700d0b770aSPeter Brune } 710d0b770aSPeter Brune } 720d0b770aSPeter Brune if (Xdot) { 730d0b770aSPeter Brune if (dm && dm != ts->dm) { 740d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 750d0b770aSPeter Brune } 760d0b770aSPeter Brune } 770d0b770aSPeter Brune PetscFunctionReturn(0); 780d0b770aSPeter Brune } 790d0b770aSPeter Brune 807445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 817445fe48SJed Brown { 827445fe48SJed Brown PetscFunctionBegin; 837445fe48SJed Brown PetscFunctionReturn(0); 847445fe48SJed Brown } 857445fe48SJed Brown 867445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 877445fe48SJed Brown { 887445fe48SJed Brown TS ts = (TS)ctx; 897445fe48SJed Brown PetscErrorCode ierr; 907445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 917445fe48SJed Brown 927445fe48SJed Brown PetscFunctionBegin; 937445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 947445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 957445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 967445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 977445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 987445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 990d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 1000d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 1017445fe48SJed Brown PetscFunctionReturn(0); 1027445fe48SJed Brown } 1037445fe48SJed Brown 104258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 105258e1594SPeter Brune { 106258e1594SPeter Brune PetscFunctionBegin; 107258e1594SPeter Brune PetscFunctionReturn(0); 108258e1594SPeter Brune } 109258e1594SPeter Brune 110258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 111258e1594SPeter Brune { 112258e1594SPeter Brune TS ts = (TS)ctx; 113258e1594SPeter Brune PetscErrorCode ierr; 114258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 115258e1594SPeter Brune 116258e1594SPeter Brune PetscFunctionBegin; 117258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 118258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 119258e1594SPeter Brune 120258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 122258e1594SPeter Brune 123258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 124258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125258e1594SPeter Brune 126258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 127258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 128258e1594SPeter Brune PetscFunctionReturn(0); 129258e1594SPeter Brune } 130258e1594SPeter Brune 131022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 132022c081aSHong Zhang { 133022c081aSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 134cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 135022c081aSHong Zhang PetscErrorCode ierr; 136022c081aSHong Zhang 137022c081aSHong Zhang PetscFunctionBegin; 138022c081aSHong Zhang if (th->endpoint) { 139022c081aSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 140022c081aSHong Zhang if (th->Theta!=1.0) { 141cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 142cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 143022c081aSHong Zhang } 144cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 145cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 146022c081aSHong Zhang } else { 147cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 148cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 149022c081aSHong Zhang } 150022c081aSHong Zhang PetscFunctionReturn(0); 151022c081aSHong Zhang } 152022c081aSHong Zhang 153b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 154b1cb36f3SHong Zhang { 155b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 156cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 157b1cb36f3SHong Zhang PetscErrorCode ierr; 158b1cb36f3SHong Zhang 159b1cb36f3SHong Zhang PetscFunctionBegin; 160b1cb36f3SHong Zhang /* backup cost integral */ 161cd4cee2dSHong Zhang ierr = VecCopy(quadts->vec_sol,th->VecCostIntegral0);CHKERRQ(ierr); 162022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 163b1cb36f3SHong Zhang PetscFunctionReturn(0); 164b1cb36f3SHong Zhang } 165b1cb36f3SHong Zhang 166b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 167b1cb36f3SHong Zhang { 168b1cb36f3SHong Zhang PetscErrorCode ierr; 169b1cb36f3SHong Zhang 170b1cb36f3SHong Zhang PetscFunctionBegin; 171022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 17224655328SShri PetscFunctionReturn(0); 17324655328SShri } 17424655328SShri 175470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x) 1761566a47fSLisandro Dalcin { 1771566a47fSLisandro Dalcin PetscInt nits,lits; 1781566a47fSLisandro Dalcin PetscErrorCode ierr; 1791566a47fSLisandro Dalcin 1801566a47fSLisandro Dalcin PetscFunctionBegin; 1811566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1821566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1831566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1841566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1851566a47fSLisandro Dalcin PetscFunctionReturn(0); 1861566a47fSLisandro Dalcin } 1871566a47fSLisandro Dalcin 188193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 189316643e7SJed Brown { 190316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1911566a47fSLisandro Dalcin PetscInt rejections = 0; 1924957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 1931566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 194051f2191SLisandro Dalcin PetscErrorCode ierr; 195316643e7SJed Brown 196316643e7SJed Brown PetscFunctionBegin; 1971566a47fSLisandro Dalcin if (!ts->steprollback) { 1982ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 1993b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 2001566a47fSLisandro Dalcin } 2011566a47fSLisandro Dalcin 2021566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 2031566a47fSLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 204715f1b00SHong Zhang 205*3e05ceb1SHong Zhang th->shift = 1/(th->Theta*ts->time_step); 2061566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 207316643e7SJed Brown 208be5899b3SLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 209fecfb714SLisandro Dalcin if (th->extrapolate && !ts->steprestart) { 210*3e05ceb1SHong Zhang ierr = VecAXPY(th->X,1/th->shift,th->Xdot);CHKERRQ(ierr); 211be5899b3SLisandro Dalcin } 212eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 213eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2141566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2151566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2161566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2171566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2181566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 219eb284becSJed Brown } 220be5899b3SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 221470880abSPatrick Sanan ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 222be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 2231566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 224fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 225051f2191SLisandro Dalcin 226051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2271566a47fSLisandro Dalcin if (th->endpoint) { 2281566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2291566a47fSLisandro Dalcin } else { 230*3e05ceb1SHong Zhang ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr); 2311566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2321566a47fSLisandro Dalcin } 2331566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2341566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2351566a47fSLisandro Dalcin if (!accept) { 2361566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 237be5899b3SLisandro Dalcin ts->time_step = next_time_step; 238be5899b3SLisandro Dalcin goto reject_step; 239051f2191SLisandro Dalcin } 2403b1890cdSShri Abhyankar 241715f1b00SHong Zhang if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 242b1cb36f3SHong Zhang th->ptime = ts->ptime; 243b1cb36f3SHong Zhang th->time_step = ts->time_step; 24417145e75SHong Zhang } 2451566a47fSLisandro Dalcin 2462b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 247cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 248051f2191SLisandro Dalcin break; 249051f2191SLisandro Dalcin 250051f2191SLisandro Dalcin reject_step: 251fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 2521566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 253051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2541566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2553b1890cdSShri Abhyankar } 2563b1890cdSShri Abhyankar } 257316643e7SJed Brown PetscFunctionReturn(0); 258316643e7SJed Brown } 259316643e7SJed Brown 2606980b705SHong Zhang /* 2616980b705SHong Zhang Use SNES to compute the Jacobian so that finite differencing could be used when TS Jacobian is not available. 2626980b705SHong Zhang */ 2636980b705SHong Zhang static PetscErrorCode KSPTSFormOperator_Private(KSP ksp,Vec x,Mat J,Mat Jpre,TS ts) 2646980b705SHong Zhang { 2656980b705SHong Zhang SNES snes = ts->snes; 2666980b705SHong Zhang MatFDColoring color; 2676980b705SHong Zhang PetscErrorCode ierr; 2686980b705SHong Zhang 2696980b705SHong Zhang PetscFunctionBegin; 2706980b705SHong Zhang /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */ 2716980b705SHong Zhang ierr = PetscObjectQuery((PetscObject)Jpre,"SNESMatFDColoring",(PetscObject*)&color);CHKERRQ(ierr); 2726980b705SHong Zhang if (color) { 2736980b705SHong Zhang Vec f; 2746980b705SHong Zhang ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr); 2756980b705SHong Zhang ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); 2766980b705SHong Zhang } 2776980b705SHong Zhang ierr = SNESComputeJacobian(snes,x,J,Jpre);CHKERRQ(ierr); 2786980b705SHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 2796980b705SHong Zhang PetscFunctionReturn(0); 2806980b705SHong Zhang } 2816980b705SHong Zhang 282c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 283c9681e5dSHong Zhang { 284c9681e5dSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 285cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 286c9681e5dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 287c9681e5dSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 288c9681e5dSHong Zhang PetscInt nadj; 2897207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 290c9681e5dSHong Zhang KSP ksp; 291c9681e5dSHong Zhang PetscScalar *xarr; 292077c4feaSHong Zhang TSEquationType eqtype; 293077c4feaSHong Zhang PetscBool isexplicitode = PETSC_FALSE; 294c9681e5dSHong Zhang PetscErrorCode ierr; 295c9681e5dSHong Zhang 296c9681e5dSHong Zhang PetscFunctionBegin; 297077c4feaSHong Zhang ierr = TSGetEquationType(ts,&eqtype);CHKERRQ(ierr); 298077c4feaSHong Zhang if (eqtype == TS_EQ_ODE_EXPLICIT) { 299077c4feaSHong Zhang isexplicitode = PETSC_TRUE; 300077c4feaSHong Zhang VecsDeltaLam = ts->vecs_sensi; 301077c4feaSHong Zhang VecsDeltaLam2 = ts->vecs_sensi2; 302077c4feaSHong Zhang } 303c9681e5dSHong Zhang th->status = TS_STEP_INCOMPLETE; 304c9681e5dSHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 305cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 3067207e0fdSHong Zhang if (quadts) { 307cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 308cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 3097207e0fdSHong Zhang } 310c9681e5dSHong Zhang 311c9681e5dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 312c9681e5dSHong Zhang th->stage_time = ts->ptime; /* time_step is negative*/ 313c9681e5dSHong Zhang th->ptime = ts->ptime + ts->time_step; 314c9681e5dSHong Zhang th->time_step = -ts->time_step; 315c9681e5dSHong Zhang 31633f72c5dSHong Zhang /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */ 3177207e0fdSHong Zhang if (quadts) { 318cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 3197207e0fdSHong Zhang } 320cd4cee2dSHong Zhang 321c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 322c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 323c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */ 324cd4cee2dSHong Zhang if (quadJ) { 325cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 326cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 327cd4cee2dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr); 328cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 329cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 330c9681e5dSHong Zhang } 331c9681e5dSHong Zhang } 332c9681e5dSHong Zhang 333c9681e5dSHong Zhang /* Build LHS for first-order adjoint */ 3346980b705SHong Zhang ierr = KSPTSFormOperator_Private(ksp,th->X,J,Jpre,ts);CHKERRQ(ierr); 335c9681e5dSHong Zhang 336c9681e5dSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 337c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 338c9681e5dSHong Zhang KSPConvergedReason kspreason; 339c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 340c9681e5dSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 341c9681e5dSHong Zhang if (kspreason < 0) { 342c9681e5dSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 343c9681e5dSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 344c9681e5dSHong Zhang } 345c9681e5dSHong Zhang } 346c9681e5dSHong Zhang 347c9681e5dSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 348c9681e5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 349c9681e5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 350c9681e5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 351c9681e5dSHong Zhang /* lambda_s^T F_UU w_1 */ 352b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 353c9681e5dSHong Zhang /* lambda_s^T F_UP w_2 */ 354b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 355c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 356c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 3576980b705SHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],1./th->time_step);CHKERRQ(ierr); 35833f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 359c9681e5dSHong Zhang if (ts->vecs_fup) { 36033f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 361c9681e5dSHong Zhang } 362c9681e5dSHong Zhang } 363c9681e5dSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 364c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 36505c9950eSHong Zhang KSPConvergedReason kspreason; 366c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 36705c9950eSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 36805c9950eSHong Zhang if (kspreason < 0) { 36905c9950eSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 37005c9950eSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 37105c9950eSHong Zhang } 372c9681e5dSHong Zhang } 373c9681e5dSHong Zhang } 374c9681e5dSHong Zhang 375c9681e5dSHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 376077c4feaSHong Zhang if (!isexplicitode) { 377*3e05ceb1SHong Zhang th->shift = 0.0; 378*3e05ceb1SHong Zhang ierr = KSPTSFormOperator_Private(ksp,th->X,J,Jpre,ts);CHKERRQ(ierr); 379c9681e5dSHong Zhang ierr = MatScale(J,-1.);CHKERRQ(ierr); 380c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 38133f72c5dSHong Zhang /* Add f_U \lambda_s to the original RHS */ 382c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 383c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr); 384c9681e5dSHong Zhang ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 385c9681e5dSHong Zhang if (ts->vecs_sensi2) { 386c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 387c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr); 388c9681e5dSHong Zhang ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 389c9681e5dSHong Zhang } 390c9681e5dSHong Zhang } 391077c4feaSHong Zhang } 392c9681e5dSHong Zhang if (ts->vecs_sensip) { 393*3e05ceb1SHong Zhang th->shift = 1./th->time_step;; 394*3e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */ 3957207e0fdSHong Zhang if (quadts) { 396cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 3977207e0fdSHong Zhang } 398c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 399c9681e5dSHong Zhang /* lambda_s^T F_PU w_1 */ 400b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 40133f72c5dSHong Zhang /* lambda_s^T F_PP w_2 */ 402b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 403c9681e5dSHong Zhang } 404cd4cee2dSHong Zhang 405c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 406c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 407c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 408cd4cee2dSHong Zhang if (quadJp) { 409cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 410cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 411cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr); 412cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 413cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 41433f72c5dSHong Zhang } 415c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 416c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 417c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 418c9681e5dSHong Zhang if (ts->vecs_fpu) { 419c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 420c9681e5dSHong Zhang } 421c9681e5dSHong Zhang if (ts->vecs_fpp) { 422c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 423c9681e5dSHong Zhang } 424c9681e5dSHong Zhang } 425c9681e5dSHong Zhang } 426c9681e5dSHong Zhang } 427c9681e5dSHong Zhang 428c9681e5dSHong Zhang if (ts->vecs_sensi2) { 429c9681e5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 430c9681e5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 431c9681e5dSHong Zhang } 432c9681e5dSHong Zhang th->status = TS_STEP_COMPLETE; 433c9681e5dSHong Zhang PetscFunctionReturn(0); 434c9681e5dSHong Zhang } 435c9681e5dSHong Zhang 43642f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 4372ca6e920SHong Zhang { 4382ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 439cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 440b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 441b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 4422ca6e920SHong Zhang PetscInt nadj; 4437207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 4442ca6e920SHong Zhang KSP ksp; 445b5b4017aSHong Zhang PetscScalar *xarr; 446b5b4017aSHong Zhang PetscErrorCode ierr; 4472ca6e920SHong Zhang 4482ca6e920SHong Zhang PetscFunctionBegin; 449c9681e5dSHong Zhang if (th->Theta == 1.) { 450c9681e5dSHong Zhang ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr); 451c9681e5dSHong Zhang PetscFunctionReturn(0); 452c9681e5dSHong Zhang } 4532ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 454302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 455cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 4567207e0fdSHong Zhang if (quadts) { 457cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 458cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 4597207e0fdSHong Zhang } 460b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 461022c081aSHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 462b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 463022c081aSHong Zhang th->time_step = -ts->time_step; 4642ca6e920SHong Zhang 465b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 46633f72c5dSHong Zhang /* Cost function has an integral term */ 4677207e0fdSHong Zhang if (quadts) { 46805755b9cSHong Zhang if (th->endpoint) { 469cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr); 47005755b9cSHong Zhang } else { 471cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 47205755b9cSHong Zhang } 4737207e0fdSHong Zhang } 47433f72c5dSHong Zhang 475abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 476b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 477022c081aSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr); 478cd4cee2dSHong Zhang if (quadJ) { 479cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 480cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 481cd4cee2dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr); 482cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 483cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 48436eaed60SHong Zhang } 4852ca6e920SHong Zhang } 4863c54f92cSHong Zhang 487b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 488*3e05ceb1SHong Zhang th->shift = 1./(th->Theta*th->time_step); 4893c54f92cSHong Zhang if (th->endpoint) { 490*3e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 4913c54f92cSHong Zhang } else { 492*3e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 4933c54f92cSHong Zhang } 494cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 4952ca6e920SHong Zhang 496b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 497abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 4981d14d03bSHong Zhang KSPConvergedReason kspreason; 499b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 5001d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 5011d14d03bSHong Zhang if (kspreason < 0) { 5021d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 5031d14d03bSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 5041d14d03bSHong Zhang } 5052ca6e920SHong Zhang } 5063c54f92cSHong Zhang 5071d14d03bSHong Zhang /* Second-order adjoint */ 508b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 509b5b4017aSHong Zhang if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 510b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 511b5b4017aSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 512b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 513b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 514b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 515b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 516b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 517b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 518b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 519b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 520b5b4017aSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 521*3e05ceb1SHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],th->shift);CHKERRQ(ierr); 52233f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 523b5b4017aSHong Zhang if (ts->vecs_fup) { 52433f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 525b5b4017aSHong Zhang } 526b5b4017aSHong Zhang } 527b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 528b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 5291d14d03bSHong Zhang KSPConvergedReason kspreason; 5301d14d03bSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 5311d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 5321d14d03bSHong Zhang if (kspreason < 0) { 5331d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 5341d14d03bSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 5351d14d03bSHong Zhang } 536b5b4017aSHong Zhang } 537b5b4017aSHong Zhang } 538b5b4017aSHong Zhang 53936eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 5401d14d03bSHong Zhang if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 541*3e05ceb1SHong Zhang th->shift = 1./((th->Theta-1.)*th->time_step); 542*3e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 54333f72c5dSHong Zhang /* R_U at t_n */ 5447207e0fdSHong Zhang if (quadts) { 545cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr); 5467207e0fdSHong Zhang } 547abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 548b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 549cd4cee2dSHong Zhang if (quadJ) { 550cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 551cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 552cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr); 553cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 554cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 555b5b4017aSHong Zhang } 556*3e05ceb1SHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./th->shift);CHKERRQ(ierr); 557b5b4017aSHong Zhang } 5581d14d03bSHong Zhang 5591d14d03bSHong Zhang /* Second-order adjoint */ 5601d14d03bSHong Zhang if (ts->vecs_sensi2) { /* U_n */ 561b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 562b5b4017aSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 563b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 564b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 565b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 566b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 56733f72c5dSHong Zhang ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr); 568b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 569b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 570b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 57133f72c5dSHong Zhang /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2 */ 572b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 57333f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 574b5b4017aSHong Zhang if (ts->vecs_fup) { 57533f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 576b5b4017aSHong Zhang } 577*3e05ceb1SHong Zhang ierr = VecScale(ts->vecs_sensi2[nadj],1./th->shift);CHKERRQ(ierr); 578b5b4017aSHong Zhang } 579b5e0532dSHong Zhang } 5803fd52205SHong Zhang 5813fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 582b5b4017aSHong Zhang /* U_{n+1} */ 583*3e05ceb1SHong Zhang th->shift = -1./(th->Theta*th->time_step); 584*3e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 5857207e0fdSHong Zhang if (quadts) { 586cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 5877207e0fdSHong Zhang } 588abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 589b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 59033f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 5913fd52205SHong Zhang } 59233f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 59333f72c5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 59433f72c5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 59533f72c5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 596b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 597b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 59833f72c5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 59933f72c5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 60033f72c5dSHong Zhang 601b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 602b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 603b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 60433f72c5dSHong Zhang /* Mu2 <- Mu2 + h theta F_P^T Lambda_s + h theta (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */ 605b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 60633f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 607b5b4017aSHong Zhang if (ts->vecs_fpu) { 60833f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 609b5b4017aSHong Zhang } 610b5b4017aSHong Zhang if (ts->vecs_fpp) { 61133f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 612b5b4017aSHong Zhang } 613b5b4017aSHong Zhang } 614b5b4017aSHong Zhang } 615b5b4017aSHong Zhang 616b5b4017aSHong Zhang /* U_s */ 617*3e05ceb1SHong Zhang th->shift = 1./((th->Theta-1.0)*th->time_step); 618*3e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 6197207e0fdSHong Zhang if (quadts) { 620cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr); 6217207e0fdSHong Zhang } 622abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 623b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 62433f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 62533f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 62633f72c5dSHong Zhang /* Get w1 at t_n from TLM matrix */ 62733f72c5dSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 62833f72c5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 629b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 630b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 63133f72c5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 63233f72c5dSHong Zhang ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr); 633b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 634b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 635b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 63633f72c5dSHong Zhang /* Mu2 <- Mu2 + h(1-theta) F_P^T Lambda_s + h(1-theta) (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */ 637b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 63833f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 639b5b4017aSHong Zhang if (ts->vecs_fpu) { 64033f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 641b5e0532dSHong Zhang } 642b5b4017aSHong Zhang if (ts->vecs_fpp) { 64333f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 644b5b4017aSHong Zhang } 64536eaed60SHong Zhang } 64636eaed60SHong Zhang } 6473fd52205SHong Zhang } 648b5e0532dSHong Zhang } 6493fd52205SHong Zhang } else { /* one-stage case */ 650*3e05ceb1SHong Zhang th->shift = 0.0; 651*3e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 6527207e0fdSHong Zhang if (quadts) { 653cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 6547207e0fdSHong Zhang } 655abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 656b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 657022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 658cd4cee2dSHong Zhang if (quadJ) { 659cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 660cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 661cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vec_drdu_col);CHKERRQ(ierr); 662cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 663cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 66436eaed60SHong Zhang } 6652ca6e920SHong Zhang } 6663fd52205SHong Zhang if (ts->vecs_sensip) { 667*3e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 6687207e0fdSHong Zhang if (quadts) { 669cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 6707207e0fdSHong Zhang } 671abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 67233f72c5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 67333f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 674cd4cee2dSHong Zhang if (quadJp) { 675cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 676cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 677cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr); 678cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 679cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 68036eaed60SHong Zhang } 68136eaed60SHong Zhang } 6823fd52205SHong Zhang } 6832ca6e920SHong Zhang } 6842ca6e920SHong Zhang 6852ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6862ca6e920SHong Zhang PetscFunctionReturn(0); 6872ca6e920SHong Zhang } 6882ca6e920SHong Zhang 689cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 690cd652676SJed Brown { 691cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 692be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 693cd652676SJed Brown PetscErrorCode ierr; 694cd652676SJed Brown 695cd652676SJed Brown PetscFunctionBegin; 696a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 697be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 698be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 699cd652676SJed Brown PetscFunctionReturn(0); 700cd652676SJed Brown } 701cd652676SJed Brown 7021566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 7031566a47fSLisandro Dalcin { 7041566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 7051566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 7061566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 7077453f775SEmil Constantinescu PetscReal wltea,wlter; 7081566a47fSLisandro Dalcin PetscErrorCode ierr; 7091566a47fSLisandro Dalcin 7101566a47fSLisandro Dalcin PetscFunctionBegin; 7112ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 7121566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 713fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 7141566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 715fecfb714SLisandro Dalcin { 716be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 717be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 7181566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 7191566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 7201566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 7211566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 7221566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 7237453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 7241566a47fSLisandro Dalcin } 7251566a47fSLisandro Dalcin if (order) *order = 2; 7261566a47fSLisandro Dalcin PetscFunctionReturn(0); 7271566a47fSLisandro Dalcin } 7281566a47fSLisandro Dalcin 7291566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 7301566a47fSLisandro Dalcin { 7311566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 7327207e0fdSHong Zhang TS quadts = ts->quadraturets; 7331566a47fSLisandro Dalcin PetscErrorCode ierr; 7341566a47fSLisandro Dalcin 7351566a47fSLisandro Dalcin PetscFunctionBegin; 7361566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 737cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 738cd4cee2dSHong Zhang ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr); 7391566a47fSLisandro Dalcin } 740715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 74113b1b0ffSHong Zhang if (ts->mat_sensip) { 74213b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7436f92202bSHong Zhang } 7447207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 7457207e0fdSHong Zhang ierr = MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7466f92202bSHong Zhang } 747715f1b00SHong Zhang PetscFunctionReturn(0); 748715f1b00SHong Zhang } 749715f1b00SHong Zhang 750715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 751715f1b00SHong Zhang { 752715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 753cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 75413b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 755b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 7567207e0fdSHong Zhang PetscInt ntlm; 757715f1b00SHong Zhang KSP ksp; 7587207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 75913b1b0ffSHong Zhang PetscScalar *barr,*xarr; 760715f1b00SHong Zhang PetscErrorCode ierr; 761715f1b00SHong Zhang 762715f1b00SHong Zhang PetscFunctionBegin; 76313b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 76413b1b0ffSHong Zhang 7657207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 7667207e0fdSHong Zhang ierr = MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7676f92202bSHong Zhang } 768715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 769cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 7707207e0fdSHong Zhang if (quadts) { 771cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 772cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 7737207e0fdSHong Zhang } 774715f1b00SHong Zhang 775715f1b00SHong Zhang /* Build RHS */ 776715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 777*3e05ceb1SHong Zhang th->shift = 1./((th->Theta-1.)*th->time_step); 778*3e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 77913b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 78013b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 781715f1b00SHong Zhang 782022c081aSHong Zhang /* Add the f_p forcing terms */ 7830e11c21fSHong Zhang if (ts->Jacp) { 784*3e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 78533f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 786*3e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 78733f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 7880e11c21fSHong Zhang } 789715f1b00SHong Zhang } else { /* 1-stage method */ 790*3e05ceb1SHong Zhang th->shift = 0.0; 791*3e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 79213b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 79313b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 79413b1b0ffSHong Zhang 795022c081aSHong Zhang /* Add the f_p forcing terms */ 7960e11c21fSHong Zhang if (ts->Jacp) { 797*3e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 79833f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 799715f1b00SHong Zhang } 8000e11c21fSHong Zhang } 801715f1b00SHong Zhang 802715f1b00SHong Zhang /* Build LHS */ 803*3e05ceb1SHong Zhang th->shift = 1/(th->Theta*th->time_step); 804715f1b00SHong Zhang if (th->endpoint) { 805*3e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 806715f1b00SHong Zhang } else { 807*3e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 808715f1b00SHong Zhang } 809cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 810715f1b00SHong Zhang 811715f1b00SHong Zhang /* 812715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 813c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 814715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 815715f1b00SHong Zhang */ 816715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 8177207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 818cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr); 819cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr); 8207207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8217207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8227207e0fdSHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 823715f1b00SHong Zhang } 824715f1b00SHong Zhang } 825715f1b00SHong Zhang 826715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 827022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 828b5b4017aSHong Zhang KSPConvergedReason kspreason; 82913b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 830b5b4017aSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 831715f1b00SHong Zhang if (th->endpoint) { 83213b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 833b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 834b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 835b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 83613b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 837715f1b00SHong Zhang } else { 838b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 839715f1b00SHong Zhang } 840b5b4017aSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 841b5b4017aSHong Zhang if (kspreason < 0) { 842b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 843b5b4017aSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 844b5b4017aSHong Zhang } 845b5b4017aSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 84613b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 847715f1b00SHong Zhang } 848715f1b00SHong Zhang 849b5b4017aSHong Zhang 85013b1b0ffSHong Zhang /* 85113b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 852c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 85313b1b0ffSHong Zhang */ 8547207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 85513b1b0ffSHong Zhang if (!th->endpoint) { 8567207e0fdSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* stage sensitivity */ 857cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 858cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 8597207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8607207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8617207e0fdSHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 86213b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 86313b1b0ffSHong Zhang } else { 864cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr); 865cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 8667207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8677207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8687207e0fdSHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 86913b1b0ffSHong Zhang } 87013b1b0ffSHong Zhang } else { 87113b1b0ffSHong Zhang if (!th->endpoint) { 87213b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 873715f1b00SHong Zhang } 874715f1b00SHong Zhang } 8751566a47fSLisandro Dalcin PetscFunctionReturn(0); 8761566a47fSLisandro Dalcin } 8771566a47fSLisandro Dalcin 8786affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 8796affb6f8SHong Zhang { 8806affb6f8SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 8816affb6f8SHong Zhang 8826affb6f8SHong Zhang PetscFunctionBegin; 8836affb6f8SHong Zhang if (ns) *ns = 1; 8846affb6f8SHong Zhang if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 8856affb6f8SHong Zhang PetscFunctionReturn(0); 8866affb6f8SHong Zhang } 8876affb6f8SHong Zhang 888316643e7SJed Brown /*------------------------------------------------------------*/ 889277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 890316643e7SJed Brown { 891316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 892316643e7SJed Brown PetscErrorCode ierr; 893316643e7SJed Brown 894316643e7SJed Brown PetscFunctionBegin; 8956bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 8966bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 8973b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 898eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 8991566a47fSLisandro Dalcin 9001566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 9011566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 9021566a47fSLisandro Dalcin 903b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 904277b19d0SLisandro Dalcin PetscFunctionReturn(0); 905277b19d0SLisandro Dalcin } 906277b19d0SLisandro Dalcin 907ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts) 908ece46509SHong Zhang { 909ece46509SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 910ece46509SHong Zhang PetscErrorCode ierr; 911ece46509SHong Zhang 912ece46509SHong Zhang PetscFunctionBegin; 913ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 914ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 915ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 916ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 917ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 918ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 919ece46509SHong Zhang PetscFunctionReturn(0); 920ece46509SHong Zhang } 921ece46509SHong Zhang 922277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 923277b19d0SLisandro Dalcin { 924277b19d0SLisandro Dalcin PetscErrorCode ierr; 925277b19d0SLisandro Dalcin 926277b19d0SLisandro Dalcin PetscFunctionBegin; 927277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 928b3a6b972SJed Brown if (ts->dm) { 929b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 930b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 931b3a6b972SJed Brown } 932277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 933bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 934bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 935bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 936bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 937316643e7SJed Brown PetscFunctionReturn(0); 938316643e7SJed Brown } 939316643e7SJed Brown 940316643e7SJed Brown /* 941316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9422b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 943316643e7SJed Brown */ 9440f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 945316643e7SJed Brown { 946316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 947316643e7SJed Brown PetscErrorCode ierr; 9487445fe48SJed Brown Vec X0,Xdot; 9497445fe48SJed Brown DM dm,dmsave; 950*3e05ceb1SHong Zhang PetscReal shift = th->shift; 951316643e7SJed Brown 952316643e7SJed Brown PetscFunctionBegin; 9537445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 9545a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9557445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 956b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 9577445fe48SJed Brown 9587445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9597445fe48SJed Brown dmsave = ts->dm; 9607445fe48SJed Brown ts->dm = dm; 9617445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 9627445fe48SJed Brown ts->dm = dmsave; 9630d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 964316643e7SJed Brown PetscFunctionReturn(0); 965316643e7SJed Brown } 966316643e7SJed Brown 967d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 968316643e7SJed Brown { 969316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 970316643e7SJed Brown PetscErrorCode ierr; 9717445fe48SJed Brown Vec Xdot; 9727445fe48SJed Brown DM dm,dmsave; 973*3e05ceb1SHong Zhang PetscReal shift = th->shift; 974316643e7SJed Brown 975316643e7SJed Brown PetscFunctionBegin; 9767445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 977be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 9780298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 9797445fe48SJed Brown 9807445fe48SJed Brown dmsave = ts->dm; 9817445fe48SJed Brown ts->dm = dm; 982d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 9837445fe48SJed Brown ts->dm = dmsave; 9840298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 985316643e7SJed Brown PetscFunctionReturn(0); 986316643e7SJed Brown } 987316643e7SJed Brown 988715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 989715f1b00SHong Zhang { 990715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 9917207e0fdSHong Zhang TS quadts = ts->quadraturets; 992715f1b00SHong Zhang PetscErrorCode ierr; 993715f1b00SHong Zhang 994715f1b00SHong Zhang PetscFunctionBegin; 995022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 99613b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 99713b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 9987207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 9997207e0fdSHong Zhang ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 10007207e0fdSHong Zhang ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr); 1001715f1b00SHong Zhang } 10026f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 100313b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 100413b1b0ffSHong Zhang 1005b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 1006715f1b00SHong Zhang PetscFunctionReturn(0); 1007715f1b00SHong Zhang } 1008715f1b00SHong Zhang 10097adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts) 10107adebddeSHong Zhang { 10117adebddeSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 10127207e0fdSHong Zhang TS quadts = ts->quadraturets; 10137adebddeSHong Zhang PetscErrorCode ierr; 10147adebddeSHong Zhang 10157adebddeSHong Zhang PetscFunctionBegin; 10167207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 10177207e0fdSHong Zhang ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr); 10187207e0fdSHong Zhang ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr); 10197adebddeSHong Zhang } 10207adebddeSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 10217adebddeSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 10227adebddeSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 10237adebddeSHong Zhang PetscFunctionReturn(0); 10247adebddeSHong Zhang } 10257adebddeSHong Zhang 1026316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 1027316643e7SJed Brown { 1028316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1029cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 10302ffb9264SLisandro Dalcin PetscBool match; 1031316643e7SJed Brown PetscErrorCode ierr; 1032316643e7SJed Brown 1033316643e7SJed Brown PetscFunctionBegin; 1034cd4cee2dSHong Zhang if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */ 1035cd4cee2dSHong Zhang ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr); 1036b1cb36f3SHong Zhang } 103739d13666SHong Zhang if (!th->X) { 1038316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 103939d13666SHong Zhang } 104039d13666SHong Zhang if (!th->Xdot) { 1041316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 104239d13666SHong Zhang } 104339d13666SHong Zhang if (!th->X0) { 10443b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 104539d13666SHong Zhang } 10461566a47fSLisandro Dalcin if (th->endpoint) { 10471566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 10487445fe48SJed Brown } 10493b1890cdSShri Abhyankar 10501566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 10511566a47fSLisandro Dalcin 10521566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 10531566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10541566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10551566a47fSLisandro Dalcin 10561566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 10571566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 10582ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 10592ffb9264SLisandro Dalcin if (!match) { 10601566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 10611566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 10623b1890cdSShri Abhyankar } 10631566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 1064b4dd3bf9SBarry Smith PetscFunctionReturn(0); 1065b4dd3bf9SBarry Smith } 10660c7376e5SHong Zhang 1067b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1068b4dd3bf9SBarry Smith 1069b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1070b4dd3bf9SBarry Smith { 1071b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 1072b4dd3bf9SBarry Smith PetscErrorCode ierr; 1073b4dd3bf9SBarry Smith 1074b4dd3bf9SBarry Smith PetscFunctionBegin; 1075b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 1076b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 10772ca6e920SHong Zhang if (ts->vecs_sensip) { 1078b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 10792ca6e920SHong Zhang } 1080b5b4017aSHong Zhang if (ts->vecs_sensi2) { 1081b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 1082b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 108367633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 108467633408SHong Zhang if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu; 108567633408SHong Zhang if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup; 1086b5b4017aSHong Zhang } 1087c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 1088c9681e5dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 108967633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 109067633408SHong Zhang if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu; 109167633408SHong Zhang if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp; 1092b5b4017aSHong Zhang } 1093316643e7SJed Brown PetscFunctionReturn(0); 1094316643e7SJed Brown } 1095316643e7SJed Brown 10964416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 1097316643e7SJed Brown { 1098316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1099316643e7SJed Brown PetscErrorCode ierr; 1100316643e7SJed Brown 1101316643e7SJed Brown PetscFunctionBegin; 1102e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 1103316643e7SJed Brown { 11040298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 11050298fd71SBarry 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); 110603842d09SLisandro 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); 1107316643e7SJed Brown } 1108316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1109316643e7SJed Brown PetscFunctionReturn(0); 1110316643e7SJed Brown } 1111316643e7SJed Brown 1112316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 1113316643e7SJed Brown { 1114316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1115ace3abfcSBarry Smith PetscBool iascii; 1116316643e7SJed Brown PetscErrorCode ierr; 1117316643e7SJed Brown 1118316643e7SJed Brown PetscFunctionBegin; 1119251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1120316643e7SJed Brown if (iascii) { 11217c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 1122ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1123316643e7SJed Brown } 1124316643e7SJed Brown PetscFunctionReturn(0); 1125316643e7SJed Brown } 1126316643e7SJed Brown 1127560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 11280de4c49aSJed Brown { 11290de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11300de4c49aSJed Brown 11310de4c49aSJed Brown PetscFunctionBegin; 11320de4c49aSJed Brown *theta = th->Theta; 11330de4c49aSJed Brown PetscFunctionReturn(0); 11340de4c49aSJed Brown } 11350de4c49aSJed Brown 1136560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 11370de4c49aSJed Brown { 11380de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11390de4c49aSJed Brown 11400de4c49aSJed Brown PetscFunctionBegin; 11417c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 11420de4c49aSJed Brown th->Theta = theta; 11431566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 11440de4c49aSJed Brown PetscFunctionReturn(0); 11450de4c49aSJed Brown } 1146eb284becSJed Brown 1147560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 114826f2ff8fSLisandro Dalcin { 114926f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 115026f2ff8fSLisandro Dalcin 115126f2ff8fSLisandro Dalcin PetscFunctionBegin; 115226f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 115326f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 115426f2ff8fSLisandro Dalcin } 115526f2ff8fSLisandro Dalcin 1156560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 1157eb284becSJed Brown { 1158eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1159eb284becSJed Brown 1160eb284becSJed Brown PetscFunctionBegin; 1161eb284becSJed Brown th->endpoint = flg; 1162eb284becSJed Brown PetscFunctionReturn(0); 1163eb284becSJed Brown } 11640de4c49aSJed Brown 1165f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1166f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 1167f9c1d6abSBarry Smith { 1168f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 1169f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 11703fd8ae06SJed Brown const PetscReal one = 1.0; 1171f9c1d6abSBarry Smith 1172f9c1d6abSBarry Smith PetscFunctionBegin; 11733fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 1174f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 1175f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 1176f9c1d6abSBarry Smith PetscFunctionReturn(0); 1177f9c1d6abSBarry Smith } 1178f9c1d6abSBarry Smith #endif 1179f9c1d6abSBarry Smith 118042682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 118142682096SHong Zhang { 118242682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 118342682096SHong Zhang 118442682096SHong Zhang PetscFunctionBegin; 11851566a47fSLisandro Dalcin if (ns) *ns = 1; 11861566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 118742682096SHong Zhang PetscFunctionReturn(0); 118842682096SHong Zhang } 1189f9c1d6abSBarry Smith 1190316643e7SJed Brown /* ------------------------------------------------------------ */ 1191316643e7SJed Brown /*MC 119296f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1193316643e7SJed Brown 1194316643e7SJed Brown Level: beginner 1195316643e7SJed Brown 11964eb428fdSBarry Smith Options Database: 1197c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1198c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 119903842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 12004eb428fdSBarry Smith 1201eb284becSJed Brown Notes: 12020c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 12030c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 12044eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 12054eb428fdSBarry Smith 1206eb284becSJed Brown This method can be applied to DAE. 1207eb284becSJed Brown 1208eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 1209eb284becSJed Brown 1210eb284becSJed Brown .vb 1211eb284becSJed Brown Theta | Theta 1212eb284becSJed Brown ------------- 1213eb284becSJed Brown | 1 1214eb284becSJed Brown .ve 1215eb284becSJed Brown 1216eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1217eb284becSJed Brown 1218eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1219eb284becSJed Brown 1220eb284becSJed Brown .vb 1221eb284becSJed Brown 0 | 0 0 1222eb284becSJed Brown 1 | 1-Theta Theta 1223eb284becSJed Brown ------------------- 1224eb284becSJed Brown | 1-Theta Theta 1225eb284becSJed Brown .ve 1226eb284becSJed Brown 1227eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1228eb284becSJed Brown 1229eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1230eb284becSJed Brown 1231eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1232eb284becSJed Brown 12334eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1234eb284becSJed Brown 1235eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1236316643e7SJed Brown 1237316643e7SJed Brown M*/ 12388cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1239316643e7SJed Brown { 1240316643e7SJed Brown TS_Theta *th; 1241316643e7SJed Brown PetscErrorCode ierr; 1242316643e7SJed Brown 1243316643e7SJed Brown PetscFunctionBegin; 1244277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1245ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1246316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1247316643e7SJed Brown ts->ops->view = TSView_Theta; 1248316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 124942f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1250ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1251316643e7SJed Brown ts->ops->step = TSStep_Theta; 1252cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12531566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 125424655328SShri ts->ops->rollback = TSRollBack_Theta; 1255316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12560f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12570f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1258f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1259f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1260f9c1d6abSBarry Smith #endif 126142682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 126242f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1263b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1264b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12652ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12666affb6f8SHong Zhang 1267715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12687adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1269715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12706affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1271316643e7SJed Brown 1272efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1273efd4aadfSBarry Smith 1274b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1275316643e7SJed Brown ts->data = (void*)th; 1276316643e7SJed Brown 1277715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1278715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1279715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1280b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1281715f1b00SHong Zhang 12826f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1283316643e7SJed Brown th->Theta = 0.5; 12841566a47fSLisandro Dalcin th->order = 2; 1285bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1286bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1287bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1288bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1289316643e7SJed Brown PetscFunctionReturn(0); 1290316643e7SJed Brown } 12910de4c49aSJed Brown 12920de4c49aSJed Brown /*@ 12930de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 12940de4c49aSJed Brown 12950de4c49aSJed Brown Not Collective 12960de4c49aSJed Brown 12970de4c49aSJed Brown Input Parameter: 12980de4c49aSJed Brown . ts - timestepping context 12990de4c49aSJed Brown 13000de4c49aSJed Brown Output Parameter: 13010de4c49aSJed Brown . theta - stage abscissa 13020de4c49aSJed Brown 13030de4c49aSJed Brown Note: 13040de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 13050de4c49aSJed Brown 13060de4c49aSJed Brown Level: Advanced 13070de4c49aSJed Brown 13080de4c49aSJed Brown .seealso: TSThetaSetTheta() 13090de4c49aSJed Brown @*/ 13107087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 13110de4c49aSJed Brown { 13124ac538c5SBarry Smith PetscErrorCode ierr; 13130de4c49aSJed Brown 13140de4c49aSJed Brown PetscFunctionBegin; 1315afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13160de4c49aSJed Brown PetscValidPointer(theta,2); 13174ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 13180de4c49aSJed Brown PetscFunctionReturn(0); 13190de4c49aSJed Brown } 13200de4c49aSJed Brown 13210de4c49aSJed Brown /*@ 13220de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 13230de4c49aSJed Brown 13240de4c49aSJed Brown Not Collective 13250de4c49aSJed Brown 13260de4c49aSJed Brown Input Parameter: 13270de4c49aSJed Brown + ts - timestepping context 13280de4c49aSJed Brown - theta - stage abscissa 13290de4c49aSJed Brown 13300de4c49aSJed Brown Options Database: 13310de4c49aSJed Brown . -ts_theta_theta <theta> 13320de4c49aSJed Brown 13330de4c49aSJed Brown Level: Intermediate 13340de4c49aSJed Brown 13350de4c49aSJed Brown .seealso: TSThetaGetTheta() 13360de4c49aSJed Brown @*/ 13377087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 13380de4c49aSJed Brown { 13394ac538c5SBarry Smith PetscErrorCode ierr; 13400de4c49aSJed Brown 13410de4c49aSJed Brown PetscFunctionBegin; 1342afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13434ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 13440de4c49aSJed Brown PetscFunctionReturn(0); 13450de4c49aSJed Brown } 1346f33bbcb6SJed Brown 134726f2ff8fSLisandro Dalcin /*@ 134826f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 134926f2ff8fSLisandro Dalcin 135026f2ff8fSLisandro Dalcin Not Collective 135126f2ff8fSLisandro Dalcin 135226f2ff8fSLisandro Dalcin Input Parameter: 135326f2ff8fSLisandro Dalcin . ts - timestepping context 135426f2ff8fSLisandro Dalcin 135526f2ff8fSLisandro Dalcin Output Parameter: 135626f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 135726f2ff8fSLisandro Dalcin 135826f2ff8fSLisandro Dalcin Level: Advanced 135926f2ff8fSLisandro Dalcin 136026f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 136126f2ff8fSLisandro Dalcin @*/ 136226f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 136326f2ff8fSLisandro Dalcin { 136426f2ff8fSLisandro Dalcin PetscErrorCode ierr; 136526f2ff8fSLisandro Dalcin 136626f2ff8fSLisandro Dalcin PetscFunctionBegin; 136726f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 136826f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1369163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 137026f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 137126f2ff8fSLisandro Dalcin } 137226f2ff8fSLisandro Dalcin 1373eb284becSJed Brown /*@ 1374eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1375eb284becSJed Brown 1376eb284becSJed Brown Not Collective 1377eb284becSJed Brown 1378eb284becSJed Brown Input Parameter: 1379eb284becSJed Brown + ts - timestepping context 1380eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1381eb284becSJed Brown 1382eb284becSJed Brown Options Database: 1383eb284becSJed Brown . -ts_theta_endpoint <flg> 1384eb284becSJed Brown 1385eb284becSJed Brown Level: Intermediate 1386eb284becSJed Brown 1387eb284becSJed Brown .seealso: TSTHETA, TSCN 1388eb284becSJed Brown @*/ 1389eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1390eb284becSJed Brown { 1391eb284becSJed Brown PetscErrorCode ierr; 1392eb284becSJed Brown 1393eb284becSJed Brown PetscFunctionBegin; 1394eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1395eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1396eb284becSJed Brown PetscFunctionReturn(0); 1397eb284becSJed Brown } 1398eb284becSJed Brown 1399f33bbcb6SJed Brown /* 1400f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1401f33bbcb6SJed Brown * The creation functions for these specializations are below. 1402f33bbcb6SJed Brown */ 1403f33bbcb6SJed Brown 14041566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 14051566a47fSLisandro Dalcin { 14061566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14071566a47fSLisandro Dalcin PetscErrorCode ierr; 14081566a47fSLisandro Dalcin 14091566a47fSLisandro Dalcin PetscFunctionBegin; 14101566a47fSLisandro 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"); 14111566a47fSLisandro 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"); 14121566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14131566a47fSLisandro Dalcin PetscFunctionReturn(0); 14141566a47fSLisandro Dalcin } 14151566a47fSLisandro Dalcin 1416f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1417f33bbcb6SJed Brown { 1418f33bbcb6SJed Brown PetscFunctionBegin; 1419f33bbcb6SJed Brown PetscFunctionReturn(0); 1420f33bbcb6SJed Brown } 1421f33bbcb6SJed Brown 1422f33bbcb6SJed Brown /*MC 1423f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1424f33bbcb6SJed Brown 1425f33bbcb6SJed Brown Level: beginner 1426f33bbcb6SJed Brown 14274eb428fdSBarry Smith Notes: 1428c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 14294eb428fdSBarry Smith 14301566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 14314eb428fdSBarry Smith 1432f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1433f33bbcb6SJed Brown 1434f33bbcb6SJed Brown M*/ 14358cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1436f33bbcb6SJed Brown { 1437f33bbcb6SJed Brown PetscErrorCode ierr; 1438f33bbcb6SJed Brown 1439f33bbcb6SJed Brown PetscFunctionBegin; 1440f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1441f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 14421566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 14430c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1444f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1445f33bbcb6SJed Brown PetscFunctionReturn(0); 1446f33bbcb6SJed Brown } 1447f33bbcb6SJed Brown 14481566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 14491566a47fSLisandro Dalcin { 14501566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14511566a47fSLisandro Dalcin PetscErrorCode ierr; 14521566a47fSLisandro Dalcin 14531566a47fSLisandro Dalcin PetscFunctionBegin; 14541566a47fSLisandro 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"); 14551566a47fSLisandro 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"); 14561566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14571566a47fSLisandro Dalcin PetscFunctionReturn(0); 14581566a47fSLisandro Dalcin } 14591566a47fSLisandro Dalcin 1460f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1461f33bbcb6SJed Brown { 1462f33bbcb6SJed Brown PetscFunctionBegin; 1463f33bbcb6SJed Brown PetscFunctionReturn(0); 1464f33bbcb6SJed Brown } 1465f33bbcb6SJed Brown 1466f33bbcb6SJed Brown /*MC 1467f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1468f33bbcb6SJed Brown 1469f33bbcb6SJed Brown Level: beginner 1470f33bbcb6SJed Brown 1471f33bbcb6SJed Brown Notes: 14727cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 14737cf5af47SJed Brown 14747cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1475f33bbcb6SJed Brown 1476f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1477f33bbcb6SJed Brown 1478f33bbcb6SJed Brown M*/ 14798cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1480f33bbcb6SJed Brown { 1481f33bbcb6SJed Brown PetscErrorCode ierr; 1482f33bbcb6SJed Brown 1483f33bbcb6SJed Brown PetscFunctionBegin; 1484f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1485f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1486eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 14870c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1488f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1489f33bbcb6SJed Brown PetscFunctionReturn(0); 1490f33bbcb6SJed Brown } 1491