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; 12cc4f23bcSHong Zhang Vec Stages[2]; /* Storage for stage solutions */ 13cc4f23bcSHong Zhang Vec X0, X, Xdot; /* Storage for u^n, u^n + dt a_{11} k_1, and time derivative u^{n+1}_t */ 141566a47fSLisandro Dalcin Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 151566a47fSLisandro Dalcin PetscReal Theta; 161a352fa8SHong Zhang PetscReal shift; /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */ 171566a47fSLisandro Dalcin PetscInt order; 181566a47fSLisandro Dalcin PetscBool endpoint; 191566a47fSLisandro Dalcin PetscBool extrapolate; 20715f1b00SHong Zhang TSStepStatus status; 211a352fa8SHong Zhang Vec VecCostIntegral0; /* Backup for roll-backs due to events, used by cost integral */ 221a352fa8SHong Zhang PetscReal ptime0; /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */ 231a352fa8SHong Zhang PetscReal time_step0; /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/ 241566a47fSLisandro Dalcin 25715f1b00SHong Zhang /* context for sensitivity analysis */ 26022c081aSHong Zhang PetscInt num_tlm; /* Total number of tangent linear equations */ 27b5e0532dSHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 28b5e0532dSHong Zhang Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 2913b1b0ffSHong Zhang Vec *VecsSensiTemp; /* Vector to be multiplied with Jacobian transpose */ 30cc4f23bcSHong Zhang Mat MatFwdStages[2]; /* TLM Stages */ 3113b1b0ffSHong Zhang Mat MatDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 32b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */ 3313b1b0ffSHong Zhang Mat MatFwdSensip0; /* backup for roll-backs due to events */ 347207e0fdSHong Zhang Mat MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 357207e0fdSHong Zhang Mat MatIntegralSensip0; /* backup for roll-backs due to events */ 36b5b4017aSHong Zhang Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */ 37b5b4017aSHong Zhang Vec *VecsDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */ 38b5b4017aSHong Zhang Vec *VecsSensi2Temp; /* Working vectors that holds the residual for the second-order adjoint */ 39b5b4017aSHong Zhang Vec *VecsAffine; /* Working vectors to store residuals */ 40715f1b00SHong Zhang /* context for error estimation */ 411566a47fSLisandro Dalcin Vec vec_sol_prev; 421566a47fSLisandro Dalcin Vec vec_lte_work; 43316643e7SJed Brown } TS_Theta; 44316643e7SJed Brown 45d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 46d71ae5a4SJacob Faibussowitsch { 477445fe48SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 487445fe48SJed Brown 497445fe48SJed Brown PetscFunctionBegin; 507445fe48SJed Brown if (X0) { 517445fe48SJed Brown if (dm && dm != ts->dm) { 529566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_X0", X0)); 537445fe48SJed Brown } else *X0 = ts->vec_sol; 547445fe48SJed Brown } 557445fe48SJed Brown if (Xdot) { 567445fe48SJed Brown if (dm && dm != ts->dm) { 579566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_Xdot", Xdot)); 587445fe48SJed Brown } else *Xdot = th->Xdot; 597445fe48SJed Brown } 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 617445fe48SJed Brown } 627445fe48SJed Brown 63d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 64d71ae5a4SJacob Faibussowitsch { 650d0b770aSPeter Brune PetscFunctionBegin; 660d0b770aSPeter Brune if (X0) { 6748a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_X0", X0)); 680d0b770aSPeter Brune } 690d0b770aSPeter Brune if (Xdot) { 7048a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_Xdot", Xdot)); 710d0b770aSPeter Brune } 723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 730d0b770aSPeter Brune } 740d0b770aSPeter Brune 75d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSTheta(DM fine, DM coarse, void *ctx) 76d71ae5a4SJacob Faibussowitsch { 777445fe48SJed Brown PetscFunctionBegin; 783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 797445fe48SJed Brown } 807445fe48SJed Brown 81d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSTheta(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 82d71ae5a4SJacob Faibussowitsch { 837445fe48SJed Brown TS ts = (TS)ctx; 847445fe48SJed Brown Vec X0, Xdot, X0_c, Xdot_c; 857445fe48SJed Brown 867445fe48SJed Brown PetscFunctionBegin; 879566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, fine, &X0, &Xdot)); 889566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 899566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, X0, X0_c)); 909566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Xdot, Xdot_c)); 919566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(X0_c, rscale, X0_c)); 929566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c)); 939566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, fine, &X0, &Xdot)); 949566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 967445fe48SJed Brown } 977445fe48SJed Brown 98d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSTheta(DM dm, DM subdm, void *ctx) 99d71ae5a4SJacob Faibussowitsch { 100258e1594SPeter Brune PetscFunctionBegin; 1013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 102258e1594SPeter Brune } 103258e1594SPeter Brune 104d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 105d71ae5a4SJacob Faibussowitsch { 106258e1594SPeter Brune TS ts = (TS)ctx; 107258e1594SPeter Brune Vec X0, Xdot, X0_sub, Xdot_sub; 108258e1594SPeter Brune 109258e1594SPeter Brune PetscFunctionBegin; 1109566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot)); 1119566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 112258e1594SPeter Brune 1139566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 1149566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 115258e1594SPeter Brune 1169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 1179566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 118258e1594SPeter Brune 1199566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 1209566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122258e1594SPeter Brune } 123258e1594SPeter Brune 124d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 125d71ae5a4SJacob Faibussowitsch { 126022c081aSHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 127cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 128022c081aSHong Zhang 129022c081aSHong Zhang PetscFunctionBegin; 130022c081aSHong Zhang if (th->endpoint) { 131022c081aSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 132022c081aSHong Zhang if (th->Theta != 1.0) { 1339566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, th->ptime0, th->X0, ts->vec_costintegrand)); 1349566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * (1.0 - th->Theta), ts->vec_costintegrand)); 135022c081aSHong Zhang } 1369566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, ts->ptime, ts->vec_sol, ts->vec_costintegrand)); 1379566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * th->Theta, ts->vec_costintegrand)); 138022c081aSHong Zhang } else { 1399566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, th->stage_time, th->X, ts->vec_costintegrand)); 1409566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, th->time_step0, ts->vec_costintegrand)); 141022c081aSHong Zhang } 1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143022c081aSHong Zhang } 144022c081aSHong Zhang 145d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 146d71ae5a4SJacob Faibussowitsch { 147b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 148cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 149b1cb36f3SHong Zhang 150b1cb36f3SHong Zhang PetscFunctionBegin; 151b1cb36f3SHong Zhang /* backup cost integral */ 1529566063dSJacob Faibussowitsch PetscCall(VecCopy(quadts->vec_sol, th->VecCostIntegral0)); 1539566063dSJacob Faibussowitsch PetscCall(TSThetaEvaluateCostIntegral(ts)); 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155b1cb36f3SHong Zhang } 156b1cb36f3SHong Zhang 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 158d71ae5a4SJacob Faibussowitsch { 1591a352fa8SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 160b1cb36f3SHong Zhang 161b1cb36f3SHong Zhang PetscFunctionBegin; 1621a352fa8SHong Zhang /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */ 1631a352fa8SHong Zhang th->ptime0 = ts->ptime + ts->time_step; 1641a352fa8SHong Zhang th->time_step0 = -ts->time_step; 1659566063dSJacob Faibussowitsch PetscCall(TSThetaEvaluateCostIntegral(ts)); 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16724655328SShri } 16824655328SShri 169d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTheta_SNESSolve(TS ts, Vec b, Vec x) 170d71ae5a4SJacob Faibussowitsch { 1711566a47fSLisandro Dalcin PetscInt nits, lits; 1721566a47fSLisandro Dalcin 1731566a47fSLisandro Dalcin PetscFunctionBegin; 1749566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, b, x)); 1759566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 1769566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 1779371c9d4SSatish Balay ts->snes_its += nits; 1789371c9d4SSatish Balay ts->ksp_its += lits; 1793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1801566a47fSLisandro Dalcin } 1811566a47fSLisandro Dalcin 182d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Theta(TS ts) 183d71ae5a4SJacob Faibussowitsch { 184316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 1851566a47fSLisandro Dalcin PetscInt rejections = 0; 1864957b756SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 1871566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 188316643e7SJed Brown 189316643e7SJed Brown PetscFunctionBegin; 1901566a47fSLisandro Dalcin if (!ts->steprollback) { 1919566063dSJacob Faibussowitsch if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev)); 1929566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X0)); 1931566a47fSLisandro Dalcin } 1941566a47fSLisandro Dalcin 1951566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 19699260352SHong Zhang while (!ts->reason && th->status != TS_STEP_COMPLETE) { 1973e05ceb1SHong Zhang th->shift = 1 / (th->Theta * ts->time_step); 1981566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta) * ts->time_step; 1999566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, th->X)); 20048a46eb9SPierre Jolivet if (th->extrapolate && !ts->steprestart) PetscCall(VecAXPY(th->X, 1 / th->shift, th->Xdot)); 201eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 2029566063dSJacob Faibussowitsch if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol, &th->affine)); 2039566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->Xdot)); 2049566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, th->X0, th->Xdot, th->affine, PETSC_FALSE)); 2059566063dSJacob Faibussowitsch PetscCall(VecScale(th->affine, (th->Theta - 1) / th->Theta)); 2061566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2079566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->affine)); 208eb284becSJed Brown } 2099566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 2109566063dSJacob Faibussowitsch PetscCall(TSTheta_SNESSolve(ts, th->affine, th->X)); 2119566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X)); 2129566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok)); 213fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 214051f2191SLisandro Dalcin 215051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2161566a47fSLisandro Dalcin if (th->endpoint) { 2179566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X, ts->vec_sol)); 2181566a47fSLisandro Dalcin } else { 2199566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); 2209566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot)); 2211566a47fSLisandro Dalcin } 2229566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 2231566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2241566a47fSLisandro Dalcin if (!accept) { 2259566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, ts->vec_sol)); 226be5899b3SLisandro Dalcin ts->time_step = next_time_step; 227be5899b3SLisandro Dalcin goto reject_step; 228051f2191SLisandro Dalcin } 2293b1890cdSShri Abhyankar 230715f1b00SHong Zhang if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 2311a352fa8SHong Zhang th->ptime0 = ts->ptime; 2321a352fa8SHong Zhang th->time_step0 = ts->time_step; 23317145e75SHong Zhang } 2342b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 235cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 236051f2191SLisandro Dalcin break; 237051f2191SLisandro Dalcin 238051f2191SLisandro Dalcin reject_step: 2399371c9d4SSatish Balay ts->reject++; 2409371c9d4SSatish Balay accept = PETSC_FALSE; 2411566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 242051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 24363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 2443b1890cdSShri Abhyankar } 2453b1890cdSShri Abhyankar } 2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 247316643e7SJed Brown } 248316643e7SJed Brown 249d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 250d71ae5a4SJacob Faibussowitsch { 251c9681e5dSHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 252cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 253c9681e5dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp; 254c9681e5dSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp; 255c9681e5dSHong Zhang PetscInt nadj; 2567207e0fdSHong Zhang Mat J, Jpre, quadJ = NULL, quadJp = NULL; 257c9681e5dSHong Zhang KSP ksp; 258c9681e5dSHong Zhang PetscScalar *xarr; 259077c4feaSHong Zhang TSEquationType eqtype; 260077c4feaSHong Zhang PetscBool isexplicitode = PETSC_FALSE; 2611a352fa8SHong Zhang PetscReal adjoint_time_step; 262c9681e5dSHong Zhang 263c9681e5dSHong Zhang PetscFunctionBegin; 2649566063dSJacob Faibussowitsch PetscCall(TSGetEquationType(ts, &eqtype)); 265077c4feaSHong Zhang if (eqtype == TS_EQ_ODE_EXPLICIT) { 266077c4feaSHong Zhang isexplicitode = PETSC_TRUE; 267077c4feaSHong Zhang VecsDeltaLam = ts->vecs_sensi; 268077c4feaSHong Zhang VecsDeltaLam2 = ts->vecs_sensi2; 269077c4feaSHong Zhang } 270c9681e5dSHong Zhang th->status = TS_STEP_INCOMPLETE; 2719566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(ts->snes, &ksp)); 2729566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 2737207e0fdSHong Zhang if (quadts) { 2749566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 2759566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 2767207e0fdSHong Zhang } 277c9681e5dSHong Zhang 2781a352fa8SHong Zhang th->stage_time = ts->ptime; 2791a352fa8SHong Zhang adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */ 280c9681e5dSHong Zhang 28133f72c5dSHong Zhang /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */ 2821baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 283cd4cee2dSHong Zhang 284c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 2859566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); 2869566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */ 287cd4cee2dSHong Zhang if (quadJ) { 2889566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 2899566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 2909566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col)); 2919566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdu_col)); 2929566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 293c9681e5dSHong Zhang } 294c9681e5dSHong Zhang } 295c9681e5dSHong Zhang 296c9681e5dSHong Zhang /* Build LHS for first-order adjoint */ 297c87ba875SIvan Yashchuk th->shift = 1. / adjoint_time_step; 2989566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 2999566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 300c9681e5dSHong Zhang 301c9681e5dSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 302c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 303c9681e5dSHong Zhang KSPConvergedReason kspreason; 3049566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj])); 3059566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 306c9681e5dSHong Zhang if (kspreason < 0) { 307c9681e5dSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 30863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj)); 309c9681e5dSHong Zhang } 310c9681e5dSHong Zhang } 311c9681e5dSHong Zhang 312c9681e5dSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 313c9681e5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 3149566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 3159566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 316c9681e5dSHong Zhang /* lambda_s^T F_UU w_1 */ 3179566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 318c9681e5dSHong Zhang /* lambda_s^T F_UP w_2 */ 3199566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 320c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */ 3219566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 3229566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step)); 3239566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj])); 32448a46eb9SPierre Jolivet if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj])); 325c9681e5dSHong Zhang } 326c9681e5dSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 327c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 32805c9950eSHong Zhang KSPConvergedReason kspreason; 3299566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj])); 3309566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 33105c9950eSHong Zhang if (kspreason < 0) { 33205c9950eSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 33363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj)); 33405c9950eSHong Zhang } 335c9681e5dSHong Zhang } 336c9681e5dSHong Zhang } 337c9681e5dSHong Zhang 338c9681e5dSHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 339077c4feaSHong Zhang if (!isexplicitode) { 3403e05ceb1SHong Zhang th->shift = 0.0; 3419566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 3429566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 343c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 34433f72c5dSHong Zhang /* Add f_U \lambda_s to the original RHS */ 3459566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensiTemp[nadj], -1.)); 3469566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj])); 3479566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step)); 3489566063dSJacob Faibussowitsch PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj])); 349c9681e5dSHong Zhang if (ts->vecs_sensi2) { 3509566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj])); 3519566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step)); 3529566063dSJacob Faibussowitsch PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj])); 353c9681e5dSHong Zhang } 354c9681e5dSHong Zhang } 355077c4feaSHong Zhang } 356c9681e5dSHong Zhang if (ts->vecs_sensip) { 3579566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol)); 3589566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */ 3591baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 360c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 361c9681e5dSHong Zhang /* lambda_s^T F_PU w_1 */ 3629566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 36333f72c5dSHong Zhang /* lambda_s^T F_PP w_2 */ 3649566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 365c9681e5dSHong Zhang } 366cd4cee2dSHong Zhang 367c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 3689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 3699566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj])); 370cd4cee2dSHong Zhang if (quadJp) { 3719566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 3729566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 3739566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col)); 3749566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdp_col)); 3759566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 37633f72c5dSHong Zhang } 377c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 3789566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 3799566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj])); 38048a46eb9SPierre Jolivet if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj])); 38148a46eb9SPierre Jolivet if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj])); 382c9681e5dSHong Zhang } 383c9681e5dSHong Zhang } 384c9681e5dSHong Zhang } 385c9681e5dSHong Zhang 386c9681e5dSHong Zhang if (ts->vecs_sensi2) { 3879566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 3889566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 389c9681e5dSHong Zhang } 390c9681e5dSHong Zhang th->status = TS_STEP_COMPLETE; 3913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 392c9681e5dSHong Zhang } 393c9681e5dSHong Zhang 394d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_Theta(TS ts) 395d71ae5a4SJacob Faibussowitsch { 3962ca6e920SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 397cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 398b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp; 399b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp; 4002ca6e920SHong Zhang PetscInt nadj; 4017207e0fdSHong Zhang Mat J, Jpre, quadJ = NULL, quadJp = NULL; 4022ca6e920SHong Zhang KSP ksp; 403b5b4017aSHong Zhang PetscScalar *xarr; 4041a352fa8SHong Zhang PetscReal adjoint_time_step; 405da81f932SPierre Jolivet PetscReal adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, usually ts->ptime is larger than adjoint_ptime) */ 4062ca6e920SHong Zhang 4072ca6e920SHong Zhang PetscFunctionBegin; 408c9681e5dSHong Zhang if (th->Theta == 1.) { 4099566063dSJacob Faibussowitsch PetscCall(TSAdjointStepBEuler_Private(ts)); 4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 411c9681e5dSHong Zhang } 4122ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 4139566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(ts->snes, &ksp)); 4149566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 4157207e0fdSHong Zhang if (quadts) { 4169566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 4179566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 4187207e0fdSHong Zhang } 419b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 4201a352fa8SHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step); 4211a352fa8SHong Zhang adjoint_ptime = ts->ptime + ts->time_step; 4221a352fa8SHong Zhang adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */ 4232ca6e920SHong Zhang 42482ad9101SHong Zhang if (!th->endpoint) { 4255072d23cSHong Zhang /* recover th->X0 using vec_sol and the stage value th->X */ 4269566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol)); 42782ad9101SHong Zhang } 42882ad9101SHong Zhang 429b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 43033f72c5dSHong Zhang /* Cost function has an integral term */ 4317207e0fdSHong Zhang if (quadts) { 43205755b9cSHong Zhang if (th->endpoint) { 4339566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 43405755b9cSHong Zhang } else { 4359566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 43605755b9cSHong Zhang } 4377207e0fdSHong Zhang } 43833f72c5dSHong Zhang 439abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 4409566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); 4419566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step))); 442cd4cee2dSHong Zhang if (quadJ) { 4439566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 4449566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 4459566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col)); 4469566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdu_col)); 4479566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 44836eaed60SHong Zhang } 4492ca6e920SHong Zhang } 4503c54f92cSHong Zhang 451b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 4521a352fa8SHong Zhang th->shift = 1. / (th->Theta * adjoint_time_step); 4533c54f92cSHong Zhang if (th->endpoint) { 4549566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 4553c54f92cSHong Zhang } else { 4569566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); 4573c54f92cSHong Zhang } 4589566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 4592ca6e920SHong Zhang 460b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 461abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 4621d14d03bSHong Zhang KSPConvergedReason kspreason; 4639566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj])); 4649566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 4651d14d03bSHong Zhang if (kspreason < 0) { 4661d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 46763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj)); 4681d14d03bSHong Zhang } 4692ca6e920SHong Zhang } 4703c54f92cSHong Zhang 4711d14d03bSHong Zhang /* Second-order adjoint */ 472b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 4733c633725SBarry Smith PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta"); 474b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 4759566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 4769566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 477b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 4789566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 4799566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 4809566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 481b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 4829566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 483b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */ 4849566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 4859566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift)); 4869566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj])); 48748a46eb9SPierre Jolivet if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj])); 488b5b4017aSHong Zhang } 489b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 490b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 4911d14d03bSHong Zhang KSPConvergedReason kspreason; 4929566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj])); 4939566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 4941d14d03bSHong Zhang if (kspreason < 0) { 4951d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 49663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj)); 4971d14d03bSHong Zhang } 498b5b4017aSHong Zhang } 499b5b4017aSHong Zhang } 500b5b4017aSHong Zhang 50136eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 5021d14d03bSHong Zhang if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 5031a352fa8SHong Zhang th->shift = 1. / ((th->Theta - 1.) * adjoint_time_step); 5041a352fa8SHong Zhang th->stage_time = adjoint_ptime; 5059566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre)); 5069566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 50733f72c5dSHong Zhang /* R_U at t_n */ 5081baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL)); 509abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 5109566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj])); 511cd4cee2dSHong Zhang if (quadJ) { 5129566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 5139566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 5149566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col)); 5159566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdu_col)); 5169566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 517b5b4017aSHong Zhang } 5189566063dSJacob Faibussowitsch PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift)); 519b5b4017aSHong Zhang } 5201d14d03bSHong Zhang 5211d14d03bSHong Zhang /* Second-order adjoint */ 5221d14d03bSHong Zhang if (ts->vecs_sensi2) { /* U_n */ 523b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 5249566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr)); 5259566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 526b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 5279566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 5289566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 5299566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr)); 530b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 5319566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 532b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 53333f72c5dSHong 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 */ 5349566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj])); 5359566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj])); 53648a46eb9SPierre Jolivet if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj])); 5379566063dSJacob Faibussowitsch PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift)); 538b5b4017aSHong Zhang } 539b5e0532dSHong Zhang } 5403fd52205SHong Zhang 541c586f2e8SHong Zhang th->stage_time = ts->ptime; /* recover the old value */ 542c586f2e8SHong Zhang 5433fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 544b5b4017aSHong Zhang /* U_{n+1} */ 54582ad9101SHong Zhang th->shift = 1.0 / (adjoint_time_step * th->Theta); 5469566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol)); 5479566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE)); 5481baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 549abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 5509566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 5519566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj])); 5523b711c3fSHong Zhang if (quadJp) { 5539566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 5549566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 5559566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col)); 5569566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdp_col)); 5579566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 5583b711c3fSHong Zhang } 5593fd52205SHong Zhang } 56033f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 56133f72c5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 5629566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 5639566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 564b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 5659566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 5669566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 5679566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 56833f72c5dSHong Zhang 569b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 5709566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 571b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 57233f72c5dSHong 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) */ 5739566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 5749566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj])); 57548a46eb9SPierre Jolivet if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj])); 57648a46eb9SPierre Jolivet if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj])); 577b5b4017aSHong Zhang } 578b5b4017aSHong Zhang } 579b5b4017aSHong Zhang 580b5b4017aSHong Zhang /* U_s */ 5819566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->Xdot)); 5829566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE)); 5831baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp)); 584abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 5859566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 5869566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj])); 5873b711c3fSHong Zhang if (quadJp) { 5889566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 5899566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 5909566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col)); 5919566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdp_col)); 5929566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 5933b711c3fSHong Zhang } 59433f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 59533f72c5dSHong Zhang /* Get w1 at t_n from TLM matrix */ 5969566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr)); 5979566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 598b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 5999566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 6009566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 6019566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr)); 602b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 6039566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 604b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 60533f72c5dSHong 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) */ 6069566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 6079566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj])); 60848a46eb9SPierre Jolivet if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj])); 60948a46eb9SPierre Jolivet if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj])); 61036eaed60SHong Zhang } 61136eaed60SHong Zhang } 6123fd52205SHong Zhang } 613b5e0532dSHong Zhang } 6143fd52205SHong Zhang } else { /* one-stage case */ 6153e05ceb1SHong Zhang th->shift = 0.0; 6169566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */ 6179566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 6181baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 619abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 6209566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj])); 6219566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj])); 622cd4cee2dSHong Zhang if (quadJ) { 6239566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 6249566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 6259566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col)); 6269566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdu_col)); 6279566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 62836eaed60SHong Zhang } 6292ca6e920SHong Zhang } 6303fd52205SHong Zhang if (ts->vecs_sensip) { 63182ad9101SHong Zhang th->shift = 1.0 / (adjoint_time_step * th->Theta); 6329566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); 6339566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 6341baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp)); 635abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 6369566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 6379566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj])); 638cd4cee2dSHong Zhang if (quadJp) { 6399566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 6409566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 6419566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col)); 6429566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdp_col)); 6439566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 64436eaed60SHong Zhang } 64536eaed60SHong Zhang } 6463fd52205SHong Zhang } 6472ca6e920SHong Zhang } 6482ca6e920SHong Zhang 6492ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6512ca6e920SHong Zhang } 6522ca6e920SHong Zhang 653d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X) 654d71ae5a4SJacob Faibussowitsch { 655cd652676SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 656be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 657cd652676SJed Brown 658cd652676SJed Brown PetscFunctionBegin; 6599566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X)); 660be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 6619566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, dt, th->Xdot, th->X)); 6623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 663cd652676SJed Brown } 664cd652676SJed Brown 665d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) 666d71ae5a4SJacob Faibussowitsch { 6671566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 6681566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 6691566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 6707453f775SEmil Constantinescu PetscReal wltea, wlter; 6711566a47fSLisandro Dalcin 6721566a47fSLisandro Dalcin PetscFunctionBegin; 6739371c9d4SSatish Balay if (!th->vec_sol_prev) { 6749371c9d4SSatish Balay *wlte = -1; 6753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6769371c9d4SSatish Balay } 6771566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 6789371c9d4SSatish Balay if (ts->steprestart) { 6799371c9d4SSatish Balay *wlte = -1; 6803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6819371c9d4SSatish Balay } 6821566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 683fecfb714SLisandro Dalcin { 684be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 685be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev / h; 6869371c9d4SSatish Balay PetscScalar scal[3]; 6879371c9d4SSatish Balay Vec vecs[3]; 6889371c9d4SSatish Balay scal[0] = +1 / a; 6899371c9d4SSatish Balay scal[1] = -1 / (a - 1); 6909371c9d4SSatish Balay scal[2] = +1 / (a * (a - 1)); 6919371c9d4SSatish Balay vecs[0] = X; 6929371c9d4SSatish Balay vecs[1] = th->X0; 6939371c9d4SSatish Balay vecs[2] = th->vec_sol_prev; 6949566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Y)); 6959566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y, 3, scal, vecs)); 6969566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter)); 6971566a47fSLisandro Dalcin } 6981566a47fSLisandro Dalcin if (order) *order = 2; 6993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7001566a47fSLisandro Dalcin } 7011566a47fSLisandro Dalcin 702d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Theta(TS ts) 703d71ae5a4SJacob Faibussowitsch { 7041566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 7057207e0fdSHong Zhang TS quadts = ts->quadraturets; 7061566a47fSLisandro Dalcin 7071566a47fSLisandro Dalcin PetscFunctionBegin; 7089566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, ts->vec_sol)); 70948a46eb9SPierre Jolivet if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol)); 710715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 7111baa6e33SBarry Smith if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN)); 71248a46eb9SPierre Jolivet if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN)); 7133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 714715f1b00SHong Zhang } 715715f1b00SHong Zhang 716d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_Theta(TS ts) 717d71ae5a4SJacob Faibussowitsch { 718715f1b00SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 719cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 72013b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 721b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 7227207e0fdSHong Zhang PetscInt ntlm; 723715f1b00SHong Zhang KSP ksp; 7247207e0fdSHong Zhang Mat J, Jpre, quadJ = NULL, quadJp = NULL; 72513b1b0ffSHong Zhang PetscScalar *barr, *xarr; 7261a352fa8SHong Zhang PetscReal previous_shift; 727715f1b00SHong Zhang 728715f1b00SHong Zhang PetscFunctionBegin; 7291a352fa8SHong Zhang previous_shift = th->shift; 7309566063dSJacob Faibussowitsch PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN)); 73113b1b0ffSHong Zhang 73248a46eb9SPierre Jolivet if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN)); 7339566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(ts->snes, &ksp)); 7349566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 7357207e0fdSHong Zhang if (quadts) { 7369566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 7379566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 7387207e0fdSHong Zhang } 739715f1b00SHong Zhang 740715f1b00SHong Zhang /* Build RHS */ 741715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 7421a352fa8SHong Zhang th->shift = 1. / ((th->Theta - 1.) * th->time_step0); 7439566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 7449566063dSJacob Faibussowitsch PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip)); 7459566063dSJacob Faibussowitsch PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta)); 746715f1b00SHong Zhang 747022c081aSHong Zhang /* Add the f_p forcing terms */ 7480e11c21fSHong Zhang if (ts->Jacp) { 7499566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->Xdot)); 7509566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 7519566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN)); 75282ad9101SHong Zhang th->shift = previous_shift; 7539566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol)); 7549566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 7559566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN)); 7560e11c21fSHong Zhang } 757715f1b00SHong Zhang } else { /* 1-stage method */ 7583e05ceb1SHong Zhang th->shift = 0.0; 7599566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 7609566063dSJacob Faibussowitsch PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip)); 7619566063dSJacob Faibussowitsch PetscCall(MatScale(MatDeltaFwdSensip, -1.)); 76213b1b0ffSHong Zhang 763022c081aSHong Zhang /* Add the f_p forcing terms */ 7640e11c21fSHong Zhang if (ts->Jacp) { 76582ad9101SHong Zhang th->shift = previous_shift; 7669566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); 7679566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 7689566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN)); 769715f1b00SHong Zhang } 7700e11c21fSHong Zhang } 771715f1b00SHong Zhang 772715f1b00SHong Zhang /* Build LHS */ 7731a352fa8SHong Zhang th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */ 774715f1b00SHong Zhang if (th->endpoint) { 7759566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 776715f1b00SHong Zhang } else { 7779566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 778715f1b00SHong Zhang } 7799566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 780715f1b00SHong Zhang 781715f1b00SHong Zhang /* 782715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 783c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 784715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 785715f1b00SHong Zhang */ 786715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 7877207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 7889566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL)); 7899566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp)); 7909566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp)); 7919566063dSJacob Faibussowitsch PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 7929566063dSJacob Faibussowitsch PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 793715f1b00SHong Zhang } 794715f1b00SHong Zhang } 795715f1b00SHong Zhang 796715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 797022c081aSHong Zhang for (ntlm = 0; ntlm < th->num_tlm; ntlm++) { 798b5b4017aSHong Zhang KSPConvergedReason kspreason; 7999566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr)); 8009566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr)); 801715f1b00SHong Zhang if (th->endpoint) { 8029566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr)); 8039566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 8049566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col)); 8059566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 8069566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 807715f1b00SHong Zhang } else { 8089566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol)); 809715f1b00SHong Zhang } 8109566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 811b5b4017aSHong Zhang if (kspreason < 0) { 812b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 81363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm)); 814b5b4017aSHong Zhang } 8159566063dSJacob Faibussowitsch PetscCall(VecResetArray(VecDeltaFwdSensipCol)); 8169566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr)); 817715f1b00SHong Zhang } 818715f1b00SHong Zhang 81913b1b0ffSHong Zhang /* 82013b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 821c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 82213b1b0ffSHong Zhang */ 8237207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 82413b1b0ffSHong Zhang if (!th->endpoint) { 8259566063dSJacob Faibussowitsch PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */ 8269566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 8279566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp)); 8289566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp)); 8299566063dSJacob Faibussowitsch PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 8309566063dSJacob Faibussowitsch PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 8319566063dSJacob Faibussowitsch PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); 83213b1b0ffSHong Zhang } else { 8339566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 8349566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 8359566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp)); 8369566063dSJacob Faibussowitsch PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 8379566063dSJacob Faibussowitsch PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 83813b1b0ffSHong Zhang } 83913b1b0ffSHong Zhang } else { 84048a46eb9SPierre Jolivet if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); 841715f1b00SHong Zhang } 8423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8431566a47fSLisandro Dalcin } 8441566a47fSLisandro Dalcin 845d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[]) 846d71ae5a4SJacob Faibussowitsch { 8476affb6f8SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 8486affb6f8SHong Zhang 8496affb6f8SHong Zhang PetscFunctionBegin; 8507409eceaSHong Zhang if (ns) { 8517409eceaSHong Zhang if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */ 8527409eceaSHong Zhang else *ns = 2; /* endpoint form */ 8537409eceaSHong Zhang } 8545072d23cSHong Zhang if (stagesensip) { 855cc4f23bcSHong Zhang if (!th->endpoint && th->Theta != 1.0) { 8567409eceaSHong Zhang th->MatFwdStages[0] = th->MatDeltaFwdSensip; 857cc4f23bcSHong Zhang } else { 858cc4f23bcSHong Zhang th->MatFwdStages[0] = th->MatFwdSensip0; 859cc4f23bcSHong Zhang th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */ 860cc4f23bcSHong Zhang } 861cc4f23bcSHong Zhang *stagesensip = th->MatFwdStages; 8625072d23cSHong Zhang } 8633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8646affb6f8SHong Zhang } 8656affb6f8SHong Zhang 866316643e7SJed Brown /*------------------------------------------------------------*/ 867d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Theta(TS ts) 868d71ae5a4SJacob Faibussowitsch { 869316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 870316643e7SJed Brown 871316643e7SJed Brown PetscFunctionBegin; 8729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X)); 8739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Xdot)); 8749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X0)); 8759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->affine)); 8761566a47fSLisandro Dalcin 8779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_sol_prev)); 8789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work)); 8791566a47fSLisandro Dalcin 8809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->VecCostIntegral0)); 8813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 882277b19d0SLisandro Dalcin } 883277b19d0SLisandro Dalcin 884d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_Theta(TS ts) 885d71ae5a4SJacob Faibussowitsch { 886ece46509SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 887ece46509SHong Zhang 888ece46509SHong Zhang PetscFunctionBegin; 8899566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam)); 8909566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu)); 8919566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2)); 8929566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2)); 8939566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp)); 8949566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp)); 8953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 896ece46509SHong Zhang } 897ece46509SHong Zhang 898d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Theta(TS ts) 899d71ae5a4SJacob Faibussowitsch { 900277b19d0SLisandro Dalcin PetscFunctionBegin; 9019566063dSJacob Faibussowitsch PetscCall(TSReset_Theta(ts)); 902b3a6b972SJed Brown if (ts->dm) { 9039566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts)); 9049566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts)); 905b3a6b972SJed Brown } 9069566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 9079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL)); 9089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL)); 9099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL)); 9109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL)); 9113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 912316643e7SJed Brown } 913316643e7SJed Brown 914316643e7SJed Brown /* 915316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9162b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 9170fd10508SMatthew Knepley 9180fd10508SMatthew Knepley Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true, 9190fd10508SMatthew Knepley otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is 9200fd10508SMatthew Knepley U = (U_{n+1} + U0)/2 921316643e7SJed Brown */ 922d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts) 923d71ae5a4SJacob Faibussowitsch { 924316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 9257445fe48SJed Brown Vec X0, Xdot; 9267445fe48SJed Brown DM dm, dmsave; 9273e05ceb1SHong Zhang PetscReal shift = th->shift; 928316643e7SJed Brown 929316643e7SJed Brown PetscFunctionBegin; 9309566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 9315a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9329566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot)); 933494ce9fcSHong Zhang if (x != X0) { 9349566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); 935494ce9fcSHong Zhang } else { 9369566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Xdot)); 937494ce9fcSHong Zhang } 9387445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9397445fe48SJed Brown dmsave = ts->dm; 9407445fe48SJed Brown ts->dm = dm; 9419566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE)); 9427445fe48SJed Brown ts->dm = dmsave; 9439566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 9443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 945316643e7SJed Brown } 946316643e7SJed Brown 947d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts) 948d71ae5a4SJacob Faibussowitsch { 949316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 9507445fe48SJed Brown Vec Xdot; 9517445fe48SJed Brown DM dm, dmsave; 9523e05ceb1SHong Zhang PetscReal shift = th->shift; 953316643e7SJed Brown 954316643e7SJed Brown PetscFunctionBegin; 9559566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 956be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 9579566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot)); 9587445fe48SJed Brown 9597445fe48SJed Brown dmsave = ts->dm; 9607445fe48SJed Brown ts->dm = dm; 9619566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE)); 9627445fe48SJed Brown ts->dm = dmsave; 9639566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 9643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 965316643e7SJed Brown } 966316643e7SJed Brown 967d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_Theta(TS ts) 968d71ae5a4SJacob Faibussowitsch { 969715f1b00SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 9707207e0fdSHong Zhang TS quadts = ts->quadraturets; 971715f1b00SHong Zhang 972715f1b00SHong Zhang PetscFunctionBegin; 973022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 97413b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 9759566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip)); 9767207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 9779566063dSJacob Faibussowitsch PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp)); 9789566063dSJacob Faibussowitsch PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0)); 979715f1b00SHong Zhang } 9806f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 9819566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0)); 98213b1b0ffSHong Zhang 9839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol)); 9843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 985715f1b00SHong Zhang } 986715f1b00SHong Zhang 987d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_Theta(TS ts) 988d71ae5a4SJacob Faibussowitsch { 9897adebddeSHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 9907207e0fdSHong Zhang TS quadts = ts->quadraturets; 9917adebddeSHong Zhang 9927adebddeSHong Zhang PetscFunctionBegin; 9937207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 9949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&th->MatIntegralSensipTemp)); 9959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&th->MatIntegralSensip0)); 9967adebddeSHong Zhang } 9979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol)); 9989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&th->MatDeltaFwdSensip)); 9999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&th->MatFwdSensip0)); 10003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10017adebddeSHong Zhang } 10027adebddeSHong Zhang 1003d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Theta(TS ts) 1004d71ae5a4SJacob Faibussowitsch { 1005316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 1006cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 10072ffb9264SLisandro Dalcin PetscBool match; 1008316643e7SJed Brown 1009316643e7SJed Brown PetscFunctionBegin; 1010cd4cee2dSHong Zhang if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */ 10119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0)); 1012b1cb36f3SHong Zhang } 101348a46eb9SPierre Jolivet if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X)); 101448a46eb9SPierre Jolivet if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot)); 101548a46eb9SPierre Jolivet if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0)); 101648a46eb9SPierre Jolivet if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine)); 10173b1890cdSShri Abhyankar 10181566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 101920d49056SMatthew G. Knepley th->shift = 1 / (th->Theta * ts->time_step); 10201566a47fSLisandro Dalcin 10219566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &ts->dm)); 10229566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts)); 10239566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts)); 10241566a47fSLisandro Dalcin 10259566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 10269566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 10279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match)); 10282ffb9264SLisandro Dalcin if (!match) { 10299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev)); 10309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work)); 10313b1890cdSShri Abhyankar } 10329566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &ts->snes)); 1033cc4f23bcSHong Zhang 1034cc4f23bcSHong Zhang ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE; 10353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1036b4dd3bf9SBarry Smith } 10370c7376e5SHong Zhang 1038b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1039b4dd3bf9SBarry Smith 1040d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1041d71ae5a4SJacob Faibussowitsch { 1042b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta *)ts->data; 1043b4dd3bf9SBarry Smith 1044b4dd3bf9SBarry Smith PetscFunctionBegin; 10459566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam)); 10469566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp)); 104748a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu)); 1048b5b4017aSHong Zhang if (ts->vecs_sensi2) { 10499566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2)); 10509566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp)); 105167633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 105267633408SHong Zhang if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu; 105367633408SHong Zhang if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup; 1054b5b4017aSHong Zhang } 1055c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 10569566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2)); 105767633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 105867633408SHong Zhang if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu; 105967633408SHong Zhang if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp; 1060b5b4017aSHong Zhang } 10613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1062316643e7SJed Brown } 1063316643e7SJed Brown 1064d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems *PetscOptionsObject) 1065d71ae5a4SJacob Faibussowitsch { 1066316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 1067316643e7SJed Brown 1068316643e7SJed Brown PetscFunctionBegin; 1069d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options"); 1070316643e7SJed Brown { 10719566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL)); 10729566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL)); 10739566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL)); 1074316643e7SJed Brown } 1075d0609cedSBarry Smith PetscOptionsHeadEnd(); 10763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1077316643e7SJed Brown } 1078316643e7SJed Brown 1079d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer) 1080d71ae5a4SJacob Faibussowitsch { 1081316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 1082ace3abfcSBarry Smith PetscBool iascii; 1083316643e7SJed Brown 1084316643e7SJed Brown PetscFunctionBegin; 10859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1086316643e7SJed Brown if (iascii) { 10879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Theta=%g\n", (double)th->Theta)); 10889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Extrapolation=%s\n", th->extrapolate ? "yes" : "no")); 1089316643e7SJed Brown } 10903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1091316643e7SJed Brown } 1092316643e7SJed Brown 1093d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta) 1094d71ae5a4SJacob Faibussowitsch { 10950de4c49aSJed Brown TS_Theta *th = (TS_Theta *)ts->data; 10960de4c49aSJed Brown 10970de4c49aSJed Brown PetscFunctionBegin; 10980de4c49aSJed Brown *theta = th->Theta; 10993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11000de4c49aSJed Brown } 11010de4c49aSJed Brown 1102d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta) 1103d71ae5a4SJacob Faibussowitsch { 11040de4c49aSJed Brown TS_Theta *th = (TS_Theta *)ts->data; 11050de4c49aSJed Brown 11060de4c49aSJed Brown PetscFunctionBegin; 1107cad9d221SBarry Smith PetscCheck(theta > 0 && theta <= 1, 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; 11103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11110de4c49aSJed Brown } 1112eb284becSJed Brown 1113d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint) 1114d71ae5a4SJacob Faibussowitsch { 111526f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 111626f2ff8fSLisandro Dalcin 111726f2ff8fSLisandro Dalcin PetscFunctionBegin; 111826f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 11193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 112026f2ff8fSLisandro Dalcin } 112126f2ff8fSLisandro Dalcin 1122d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg) 1123d71ae5a4SJacob Faibussowitsch { 1124eb284becSJed Brown TS_Theta *th = (TS_Theta *)ts->data; 1125eb284becSJed Brown 1126eb284becSJed Brown PetscFunctionBegin; 1127eb284becSJed Brown th->endpoint = flg; 11283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1129eb284becSJed Brown } 11300de4c49aSJed Brown 1131f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1132d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) 1133d71ae5a4SJacob Faibussowitsch { 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); 11423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1143f9c1d6abSBarry Smith } 1144f9c1d6abSBarry Smith #endif 1145f9c1d6abSBarry Smith 1146d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[]) 1147d71ae5a4SJacob Faibussowitsch { 114842682096SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 114942682096SHong Zhang 115042682096SHong Zhang PetscFunctionBegin; 11517409eceaSHong Zhang if (ns) { 11527409eceaSHong Zhang if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */ 11537409eceaSHong Zhang else *ns = 2; /* endpoint form */ 11547409eceaSHong Zhang } 11555072d23cSHong Zhang if (Y) { 1156cc4f23bcSHong Zhang if (!th->endpoint && th->Theta != 1.0) { 11577409eceaSHong Zhang th->Stages[0] = th->X; 1158cc4f23bcSHong Zhang } else { 1159cc4f23bcSHong Zhang th->Stages[0] = th->X0; 1160cc4f23bcSHong Zhang th->Stages[1] = ts->vec_sol; /* stiffly accurate */ 1161cc4f23bcSHong Zhang } 1162cc4f23bcSHong Zhang *Y = th->Stages; 11635072d23cSHong Zhang } 11643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116542682096SHong Zhang } 1166f9c1d6abSBarry Smith 1167316643e7SJed Brown /* ------------------------------------------------------------ */ 1168316643e7SJed Brown /*MC 116996f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1170316643e7SJed Brown 1171316643e7SJed Brown Level: beginner 1172316643e7SJed Brown 1173bcf0153eSBarry Smith Options Database Keys: 1174c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1175c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 117603842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 11774eb428fdSBarry Smith 1178eb284becSJed Brown Notes: 117937fdd005SBarry Smith .vb 118037fdd005SBarry Smith -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 118137fdd005SBarry Smith -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 118237fdd005SBarry Smith -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 118337fdd005SBarry Smith .ve 11844eb428fdSBarry Smith 11857409eceaSHong Zhang The endpoint variant of the Theta method and backward Euler can be applied to DAE. The midpoint variant is not suitable for DAEs because it is not stiffly accurate. 1186eb284becSJed Brown 11877409eceaSHong Zhang The midpoint variant is cast as a 1-stage implicit Runge-Kutta method. 1188eb284becSJed Brown 1189eb284becSJed Brown .vb 1190eb284becSJed Brown Theta | Theta 1191eb284becSJed Brown ------------- 1192eb284becSJed Brown | 1 1193eb284becSJed Brown .ve 1194eb284becSJed Brown 1195eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1196eb284becSJed Brown 1197eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1198eb284becSJed Brown 1199eb284becSJed Brown .vb 1200eb284becSJed Brown 0 | 0 0 1201eb284becSJed Brown 1 | 1-Theta Theta 1202eb284becSJed Brown ------------------- 1203eb284becSJed Brown | 1-Theta Theta 1204eb284becSJed Brown .ve 1205eb284becSJed Brown 1206eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1207eb284becSJed Brown 1208eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1209eb284becSJed Brown 1210eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1211eb284becSJed Brown 12124eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1213eb284becSJed Brown 1214*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()` 1215316643e7SJed Brown M*/ 1216d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1217d71ae5a4SJacob Faibussowitsch { 1218316643e7SJed Brown TS_Theta *th; 1219316643e7SJed Brown 1220316643e7SJed Brown PetscFunctionBegin; 1221277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1222ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1223316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1224316643e7SJed Brown ts->ops->view = TSView_Theta; 1225316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 122642f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1227ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1228316643e7SJed Brown ts->ops->step = TSStep_Theta; 1229cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12301566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 123124655328SShri ts->ops->rollback = TSRollBack_Theta; 1232316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12330f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12340f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1235f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1236f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1237f9c1d6abSBarry Smith #endif 123842682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 123942f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1240b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1241b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12422ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12436affb6f8SHong Zhang 1244715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12457adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1246715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12476affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1248316643e7SJed Brown 1249efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1250efd4aadfSBarry Smith 12514dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&th)); 1252316643e7SJed Brown ts->data = (void *)th; 1253316643e7SJed Brown 1254715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1255715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1256715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1257b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1258715f1b00SHong Zhang 12596f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1260316643e7SJed Brown th->Theta = 0.5; 12611566a47fSLisandro Dalcin th->order = 2; 12629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta)); 12639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta)); 12649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta)); 12659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta)); 12663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1267316643e7SJed Brown } 12680de4c49aSJed Brown 12690de4c49aSJed Brown /*@ 1270bcf0153eSBarry Smith TSThetaGetTheta - Get the abscissa of the stage in (0,1] for `TSTHETA` 12710de4c49aSJed Brown 12720de4c49aSJed Brown Not Collective 12730de4c49aSJed Brown 12740de4c49aSJed Brown Input Parameter: 12750de4c49aSJed Brown . ts - timestepping context 12760de4c49aSJed Brown 12770de4c49aSJed Brown Output Parameter: 12780de4c49aSJed Brown . theta - stage abscissa 12790de4c49aSJed Brown 1280bcf0153eSBarry Smith Level: advanced 1281bcf0153eSBarry Smith 12820de4c49aSJed Brown Note: 1283bcf0153eSBarry Smith Use of this function is normally only required to hack `TSTHETA` to use a modified integration scheme. 12840de4c49aSJed Brown 1285*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetTheta()`, `TSTHETA` 12860de4c49aSJed Brown @*/ 1287d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta) 1288d71ae5a4SJacob Faibussowitsch { 12890de4c49aSJed Brown PetscFunctionBegin; 1290afb20b64SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1291dadcf809SJacob Faibussowitsch PetscValidRealPointer(theta, 2); 1292cac4c232SBarry Smith PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta)); 12933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12940de4c49aSJed Brown } 12950de4c49aSJed Brown 12960de4c49aSJed Brown /*@ 1297bcf0153eSBarry Smith TSThetaSetTheta - Set the abscissa of the stage in (0,1] for `TSTHETA` 12980de4c49aSJed Brown 12990de4c49aSJed Brown Not Collective 13000de4c49aSJed Brown 1301d8d19677SJose E. Roman Input Parameters: 13020de4c49aSJed Brown + ts - timestepping context 13030de4c49aSJed Brown - theta - stage abscissa 13040de4c49aSJed Brown 1305bcf0153eSBarry Smith Options Database Key: 130667b8a455SSatish Balay . -ts_theta_theta <theta> - set theta 13070de4c49aSJed Brown 1308bcf0153eSBarry Smith Level: intermediate 13090de4c49aSJed Brown 1310*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaGetTheta()`, `TSTHETA`, `TSCN` 13110de4c49aSJed Brown @*/ 1312d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta) 1313d71ae5a4SJacob Faibussowitsch { 13140de4c49aSJed Brown PetscFunctionBegin; 1315afb20b64SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1316cac4c232SBarry Smith PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta)); 13173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13180de4c49aSJed Brown } 1319f33bbcb6SJed Brown 132026f2ff8fSLisandro Dalcin /*@ 1321bcf0153eSBarry Smith TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA` 132226f2ff8fSLisandro Dalcin 132326f2ff8fSLisandro Dalcin Not Collective 132426f2ff8fSLisandro Dalcin 132526f2ff8fSLisandro Dalcin Input Parameter: 132626f2ff8fSLisandro Dalcin . ts - timestepping context 132726f2ff8fSLisandro Dalcin 132826f2ff8fSLisandro Dalcin Output Parameter: 1329bcf0153eSBarry Smith . endpoint - `PETSC_TRUE` when using the endpoint variant 133026f2ff8fSLisandro Dalcin 1331bcf0153eSBarry Smith Level: advanced 133226f2ff8fSLisandro Dalcin 1333*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN` 133426f2ff8fSLisandro Dalcin @*/ 1335d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint) 1336d71ae5a4SJacob Faibussowitsch { 133726f2ff8fSLisandro Dalcin PetscFunctionBegin; 133826f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1339dadcf809SJacob Faibussowitsch PetscValidBoolPointer(endpoint, 2); 1340cac4c232SBarry Smith PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint)); 13413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 134226f2ff8fSLisandro Dalcin } 134326f2ff8fSLisandro Dalcin 1344eb284becSJed Brown /*@ 1345bcf0153eSBarry Smith TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA` 1346eb284becSJed Brown 1347eb284becSJed Brown Not Collective 1348eb284becSJed Brown 1349d8d19677SJose E. Roman Input Parameters: 1350eb284becSJed Brown + ts - timestepping context 1351bcf0153eSBarry Smith - flg - `PETSC_TRUE` to use the endpoint variant 1352eb284becSJed Brown 1353bcf0153eSBarry Smith Options Database Key: 135467b8a455SSatish Balay . -ts_theta_endpoint <flg> - use the endpoint variant 1355eb284becSJed Brown 1356bcf0153eSBarry Smith Level: intermediate 1357eb284becSJed Brown 1358*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSTHETA`, `TSCN` 1359eb284becSJed Brown @*/ 1360d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg) 1361d71ae5a4SJacob Faibussowitsch { 1362eb284becSJed Brown PetscFunctionBegin; 1363eb284becSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1364cac4c232SBarry Smith PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg)); 13653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1366eb284becSJed Brown } 1367eb284becSJed Brown 1368f33bbcb6SJed Brown /* 1369f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1370f33bbcb6SJed Brown * The creation functions for these specializations are below. 1371f33bbcb6SJed Brown */ 1372f33bbcb6SJed Brown 1373d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BEuler(TS ts) 1374d71ae5a4SJacob Faibussowitsch { 13751566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 13761566a47fSLisandro Dalcin 13771566a47fSLisandro Dalcin PetscFunctionBegin; 137808401ef6SPierre Jolivet PetscCheck(th->Theta == 1.0, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change the default value (1) of theta when using backward Euler"); 137928b400f6SJacob Faibussowitsch PetscCheck(!th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change to the endpoint form of the Theta methods when using backward Euler"); 13809566063dSJacob Faibussowitsch PetscCall(TSSetUp_Theta(ts)); 13813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13821566a47fSLisandro Dalcin } 13831566a47fSLisandro Dalcin 1384d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer) 1385d71ae5a4SJacob Faibussowitsch { 1386f33bbcb6SJed Brown PetscFunctionBegin; 13873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1388f33bbcb6SJed Brown } 1389f33bbcb6SJed Brown 1390f33bbcb6SJed Brown /*MC 1391f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1392f33bbcb6SJed Brown 1393f33bbcb6SJed Brown Level: beginner 1394f33bbcb6SJed Brown 1395bcf0153eSBarry Smith Note: 139637fdd005SBarry Smith `TSBEULER` is equivalent to `TSTHETA` with Theta=1.0 or `-ts_type theta -ts_theta_theta 1.0` 13974eb428fdSBarry Smith 1398*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA` 1399f33bbcb6SJed Brown M*/ 1400d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1401d71ae5a4SJacob Faibussowitsch { 1402f33bbcb6SJed Brown PetscFunctionBegin; 14039566063dSJacob Faibussowitsch PetscCall(TSCreate_Theta(ts)); 14049566063dSJacob Faibussowitsch PetscCall(TSThetaSetTheta(ts, 1.0)); 14059566063dSJacob Faibussowitsch PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE)); 14060c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1407f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 14083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1409f33bbcb6SJed Brown } 1410f33bbcb6SJed Brown 1411d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_CN(TS ts) 1412d71ae5a4SJacob Faibussowitsch { 14131566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 14141566a47fSLisandro Dalcin 14151566a47fSLisandro Dalcin PetscFunctionBegin; 141608401ef6SPierre Jolivet PetscCheck(th->Theta == 0.5, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change the default value (0.5) of theta when using Crank-Nicolson"); 14173c633725SBarry Smith PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change to the midpoint form of the Theta methods when using Crank-Nicolson"); 14189566063dSJacob Faibussowitsch PetscCall(TSSetUp_Theta(ts)); 14193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14201566a47fSLisandro Dalcin } 14211566a47fSLisandro Dalcin 1422d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer) 1423d71ae5a4SJacob Faibussowitsch { 1424f33bbcb6SJed Brown PetscFunctionBegin; 14253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1426f33bbcb6SJed Brown } 1427f33bbcb6SJed Brown 1428f33bbcb6SJed Brown /*MC 1429f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1430f33bbcb6SJed Brown 1431f33bbcb6SJed Brown Level: beginner 1432f33bbcb6SJed Brown 1433f33bbcb6SJed Brown Notes: 1434bcf0153eSBarry Smith `TSCN` is equivalent to `TSTHETA` with Theta=0.5 and the "endpoint" option set. I.e. 1435bcf0153eSBarry Smith .vb 1436bcf0153eSBarry Smith -ts_type theta 1437bcf0153eSBarry Smith -ts_theta_theta 0.5 1438bcf0153eSBarry Smith -ts_theta_endpoint 1439bcf0153eSBarry Smith .ve 14407cf5af47SJed Brown 1441*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`, `TSType`, 1442f33bbcb6SJed Brown M*/ 1443d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1444d71ae5a4SJacob Faibussowitsch { 1445f33bbcb6SJed Brown PetscFunctionBegin; 14469566063dSJacob Faibussowitsch PetscCall(TSCreate_Theta(ts)); 14479566063dSJacob Faibussowitsch PetscCall(TSThetaSetTheta(ts, 0.5)); 14489566063dSJacob Faibussowitsch PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE)); 14490c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1450f33bbcb6SJed Brown ts->ops->view = TSView_CN; 14513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1452f33bbcb6SJed Brown } 1453