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) { 51ac530a7eSPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_X0", X0)); 52ac530a7eSPierre Jolivet else *X0 = ts->vec_sol; 537445fe48SJed Brown } 547445fe48SJed Brown if (Xdot) { 55ac530a7eSPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_Xdot", Xdot)); 56ac530a7eSPierre Jolivet else *Xdot = th->Xdot; 577445fe48SJed Brown } 583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 597445fe48SJed Brown } 607445fe48SJed Brown 61d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 62d71ae5a4SJacob Faibussowitsch { 630d0b770aSPeter Brune PetscFunctionBegin; 640d0b770aSPeter Brune if (X0) { 6548a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_X0", X0)); 660d0b770aSPeter Brune } 670d0b770aSPeter Brune if (Xdot) { 6848a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_Xdot", Xdot)); 690d0b770aSPeter Brune } 703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 710d0b770aSPeter Brune } 720d0b770aSPeter Brune 73*2a8381b2SBarry Smith static PetscErrorCode DMCoarsenHook_TSTheta(DM fine, DM coarse, PetscCtx ctx) 74d71ae5a4SJacob Faibussowitsch { 757445fe48SJed Brown PetscFunctionBegin; 763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 777445fe48SJed Brown } 787445fe48SJed Brown 79*2a8381b2SBarry Smith static PetscErrorCode DMRestrictHook_TSTheta(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx) 80d71ae5a4SJacob Faibussowitsch { 817445fe48SJed Brown TS ts = (TS)ctx; 827445fe48SJed Brown Vec X0, Xdot, X0_c, Xdot_c; 837445fe48SJed Brown 847445fe48SJed Brown PetscFunctionBegin; 859566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, fine, &X0, &Xdot)); 869566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 879566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, X0, X0_c)); 889566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Xdot, Xdot_c)); 899566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(X0_c, rscale, X0_c)); 909566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c)); 919566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, fine, &X0, &Xdot)); 929566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 947445fe48SJed Brown } 957445fe48SJed Brown 96*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainHook_TSTheta(DM dm, DM subdm, PetscCtx ctx) 97d71ae5a4SJacob Faibussowitsch { 98258e1594SPeter Brune PetscFunctionBegin; 993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 100258e1594SPeter Brune } 101258e1594SPeter Brune 102*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx ctx) 103d71ae5a4SJacob Faibussowitsch { 104258e1594SPeter Brune TS ts = (TS)ctx; 105258e1594SPeter Brune Vec X0, Xdot, X0_sub, Xdot_sub; 106258e1594SPeter Brune 107258e1594SPeter Brune PetscFunctionBegin; 1089566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot)); 1099566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 110258e1594SPeter Brune 1119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 1129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 113258e1594SPeter Brune 1149566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 1159566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 116258e1594SPeter Brune 1179566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 1189566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120258e1594SPeter Brune } 121258e1594SPeter Brune 122d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 123d71ae5a4SJacob Faibussowitsch { 124022c081aSHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 125cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 126022c081aSHong Zhang 127022c081aSHong Zhang PetscFunctionBegin; 128022c081aSHong Zhang if (th->endpoint) { 129022c081aSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 130022c081aSHong Zhang if (th->Theta != 1.0) { 1319566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, th->ptime0, th->X0, ts->vec_costintegrand)); 1329566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * (1.0 - th->Theta), ts->vec_costintegrand)); 133022c081aSHong Zhang } 1349566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, ts->ptime, ts->vec_sol, ts->vec_costintegrand)); 1359566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * th->Theta, ts->vec_costintegrand)); 136022c081aSHong Zhang } else { 1379566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, th->stage_time, th->X, ts->vec_costintegrand)); 1389566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, th->time_step0, ts->vec_costintegrand)); 139022c081aSHong Zhang } 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141022c081aSHong Zhang } 142022c081aSHong Zhang 143d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 144d71ae5a4SJacob Faibussowitsch { 145b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 146cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 147b1cb36f3SHong Zhang 148b1cb36f3SHong Zhang PetscFunctionBegin; 149b1cb36f3SHong Zhang /* backup cost integral */ 1509566063dSJacob Faibussowitsch PetscCall(VecCopy(quadts->vec_sol, th->VecCostIntegral0)); 1519566063dSJacob Faibussowitsch PetscCall(TSThetaEvaluateCostIntegral(ts)); 1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 153b1cb36f3SHong Zhang } 154b1cb36f3SHong Zhang 155d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 156d71ae5a4SJacob Faibussowitsch { 1571a352fa8SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 158b1cb36f3SHong Zhang 159b1cb36f3SHong Zhang PetscFunctionBegin; 1601a352fa8SHong Zhang /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */ 1611a352fa8SHong Zhang th->ptime0 = ts->ptime + ts->time_step; 1621a352fa8SHong Zhang th->time_step0 = -ts->time_step; 1639566063dSJacob Faibussowitsch PetscCall(TSThetaEvaluateCostIntegral(ts)); 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16524655328SShri } 16624655328SShri 167d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTheta_SNESSolve(TS ts, Vec b, Vec x) 168d71ae5a4SJacob Faibussowitsch { 1691566a47fSLisandro Dalcin PetscInt nits, lits; 1701566a47fSLisandro Dalcin 1711566a47fSLisandro Dalcin PetscFunctionBegin; 1729566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, b, x)); 1739566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 1749566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 1759371c9d4SSatish Balay ts->snes_its += nits; 1769371c9d4SSatish Balay ts->ksp_its += lits; 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1781566a47fSLisandro Dalcin } 1791566a47fSLisandro Dalcin 18026b5f05dSStefano Zampini static PetscErrorCode TSResizeRegister_Theta(TS ts, PetscBool reg) 18126b5f05dSStefano Zampini { 18226b5f05dSStefano Zampini TS_Theta *th = (TS_Theta *)ts->data; 18326b5f05dSStefano Zampini 18426b5f05dSStefano Zampini PetscFunctionBegin; 185ecc87898SStefano Zampini if (reg) { 186ecc87898SStefano Zampini PetscCall(TSResizeRegisterVec(ts, "ts:theta:sol_prev", th->vec_sol_prev)); 187ecc87898SStefano Zampini PetscCall(TSResizeRegisterVec(ts, "ts:theta:X0", th->X0)); 188ecc87898SStefano Zampini } else { 189ecc87898SStefano Zampini PetscCall(TSResizeRetrieveVec(ts, "ts:theta:sol_prev", &th->vec_sol_prev)); 190ecc87898SStefano Zampini PetscCall(PetscObjectReference((PetscObject)th->vec_sol_prev)); 191ecc87898SStefano Zampini PetscCall(TSResizeRetrieveVec(ts, "ts:theta:X0", &th->X0)); 19226b5f05dSStefano Zampini PetscCall(PetscObjectReference((PetscObject)th->X0)); 19326b5f05dSStefano Zampini } 19426b5f05dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 19526b5f05dSStefano Zampini } 19626b5f05dSStefano Zampini 197d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Theta(TS ts) 198d71ae5a4SJacob Faibussowitsch { 199316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 2001566a47fSLisandro Dalcin PetscInt rejections = 0; 2014957b756SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 2021566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 203316643e7SJed Brown 204316643e7SJed Brown PetscFunctionBegin; 2051566a47fSLisandro Dalcin if (!ts->steprollback) { 2069566063dSJacob Faibussowitsch if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev)); 2079566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X0)); 2081566a47fSLisandro Dalcin } 2091566a47fSLisandro Dalcin 2101566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 21199260352SHong Zhang while (!ts->reason && th->status != TS_STEP_COMPLETE) { 2123e05ceb1SHong Zhang th->shift = 1 / (th->Theta * ts->time_step); 2131566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta) * ts->time_step; 2149566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, th->X)); 21548a46eb9SPierre Jolivet if (th->extrapolate && !ts->steprestart) PetscCall(VecAXPY(th->X, 1 / th->shift, th->Xdot)); 216eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 2179566063dSJacob Faibussowitsch if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol, &th->affine)); 2189566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->Xdot)); 2199566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, th->X0, th->Xdot, th->affine, PETSC_FALSE)); 2209566063dSJacob Faibussowitsch PetscCall(VecScale(th->affine, (th->Theta - 1) / th->Theta)); 221eb284becSJed Brown } 2229566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 22382d7ce88SStefano Zampini PetscCall(TSTheta_SNESSolve(ts, th->endpoint ? th->affine : NULL, th->X)); 2249566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X)); 2259566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok)); 226fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 227051f2191SLisandro Dalcin 228b7795fd7SStefano Zampini if (th->endpoint) { 2299566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X, ts->vec_sol)); 2301566a47fSLisandro Dalcin } else { 231b7795fd7SStefano Zampini PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); /* th->Xdot is needed by TSInterpolate_Theta */ 232b7795fd7SStefano Zampini if (th->Theta == 1.0) PetscCall(VecCopy(th->X, ts->vec_sol)); /* BEULER, stage already checked */ 233b7795fd7SStefano Zampini else { 2349566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot)); 235b7795fd7SStefano Zampini PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, ts->vec_sol, &stageok)); 236b7795fd7SStefano Zampini if (!stageok) { 237b7795fd7SStefano Zampini PetscCall(VecCopy(th->X0, ts->vec_sol)); 238b7795fd7SStefano Zampini goto reject_step; 2391566a47fSLisandro Dalcin } 240b7795fd7SStefano Zampini } 241b7795fd7SStefano Zampini } 242b7795fd7SStefano Zampini 243b7795fd7SStefano Zampini th->status = TS_STEP_PENDING; 2449566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 2451566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2461566a47fSLisandro Dalcin if (!accept) { 2479566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, ts->vec_sol)); 248be5899b3SLisandro Dalcin ts->time_step = next_time_step; 249be5899b3SLisandro Dalcin goto reject_step; 250051f2191SLisandro Dalcin } 2513b1890cdSShri Abhyankar 252715f1b00SHong Zhang if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 2531a352fa8SHong Zhang th->ptime0 = ts->ptime; 2541a352fa8SHong Zhang th->time_step0 = ts->time_step; 25517145e75SHong Zhang } 2562b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 257cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 258051f2191SLisandro Dalcin break; 259051f2191SLisandro Dalcin 260051f2191SLisandro Dalcin reject_step: 2619371c9d4SSatish Balay ts->reject++; 2629371c9d4SSatish Balay accept = PETSC_FALSE; 2631566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 264051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 26563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 2663b1890cdSShri Abhyankar } 2673b1890cdSShri Abhyankar } 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269316643e7SJed Brown } 270316643e7SJed Brown 271d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 272d71ae5a4SJacob Faibussowitsch { 273c9681e5dSHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 274cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 275c9681e5dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp; 276c9681e5dSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp; 277c9681e5dSHong Zhang PetscInt nadj; 2787207e0fdSHong Zhang Mat J, Jpre, quadJ = NULL, quadJp = NULL; 279c9681e5dSHong Zhang KSP ksp; 280c9681e5dSHong Zhang PetscScalar *xarr; 281077c4feaSHong Zhang TSEquationType eqtype; 282077c4feaSHong Zhang PetscBool isexplicitode = PETSC_FALSE; 2831a352fa8SHong Zhang PetscReal adjoint_time_step; 284c9681e5dSHong Zhang 285c9681e5dSHong Zhang PetscFunctionBegin; 2869566063dSJacob Faibussowitsch PetscCall(TSGetEquationType(ts, &eqtype)); 287077c4feaSHong Zhang if (eqtype == TS_EQ_ODE_EXPLICIT) { 288077c4feaSHong Zhang isexplicitode = PETSC_TRUE; 289077c4feaSHong Zhang VecsDeltaLam = ts->vecs_sensi; 290077c4feaSHong Zhang VecsDeltaLam2 = ts->vecs_sensi2; 291077c4feaSHong Zhang } 292c9681e5dSHong Zhang th->status = TS_STEP_INCOMPLETE; 2939566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(ts->snes, &ksp)); 2949566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 2957207e0fdSHong Zhang if (quadts) { 2969566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 2979566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 2987207e0fdSHong Zhang } 299c9681e5dSHong Zhang 3001a352fa8SHong Zhang th->stage_time = ts->ptime; 3011a352fa8SHong Zhang adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */ 302c9681e5dSHong Zhang 30333f72c5dSHong Zhang /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */ 3041baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 305cd4cee2dSHong Zhang 306c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 3079566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); 3089566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */ 309cd4cee2dSHong Zhang if (quadJ) { 3109566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 3119566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 3129566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col)); 3139566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdu_col)); 3149566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 315c9681e5dSHong Zhang } 316c9681e5dSHong Zhang } 317c9681e5dSHong Zhang 318c9681e5dSHong Zhang /* Build LHS for first-order adjoint */ 319c87ba875SIvan Yashchuk th->shift = 1. / adjoint_time_step; 3209566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 3219566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 322c9681e5dSHong Zhang 323c9681e5dSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 324c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 325c9681e5dSHong Zhang KSPConvergedReason kspreason; 3269566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj])); 3279566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 328c9681e5dSHong Zhang if (kspreason < 0) { 329c9681e5dSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 33063a3b9bcSJacob 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)); 331c9681e5dSHong Zhang } 332c9681e5dSHong Zhang } 333c9681e5dSHong Zhang 334c9681e5dSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 335c9681e5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 3369566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 3379566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 338c9681e5dSHong Zhang /* lambda_s^T F_UU w_1 */ 3399566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 340c9681e5dSHong Zhang /* lambda_s^T F_UP w_2 */ 3419566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 342c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */ 3439566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 3449566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step)); 3459566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj])); 34648a46eb9SPierre Jolivet if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj])); 347c9681e5dSHong Zhang } 348c9681e5dSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 349c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 35005c9950eSHong Zhang KSPConvergedReason kspreason; 3519566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj])); 3529566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 35305c9950eSHong Zhang if (kspreason < 0) { 35405c9950eSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 35563a3b9bcSJacob 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)); 35605c9950eSHong Zhang } 357c9681e5dSHong Zhang } 358c9681e5dSHong Zhang } 359c9681e5dSHong Zhang 360c9681e5dSHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 361077c4feaSHong Zhang if (!isexplicitode) { 3623e05ceb1SHong Zhang th->shift = 0.0; 3639566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 3649566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 365c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 36633f72c5dSHong Zhang /* Add f_U \lambda_s to the original RHS */ 3679566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensiTemp[nadj], -1.)); 3689566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj])); 3699566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step)); 3709566063dSJacob Faibussowitsch PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj])); 371c9681e5dSHong Zhang if (ts->vecs_sensi2) { 3729566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj])); 3739566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step)); 3749566063dSJacob Faibussowitsch PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj])); 375c9681e5dSHong Zhang } 376c9681e5dSHong Zhang } 377077c4feaSHong Zhang } 378c9681e5dSHong Zhang if (ts->vecs_sensip) { 3799566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol)); 3809566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */ 3811baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 382c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 383c9681e5dSHong Zhang /* lambda_s^T F_PU w_1 */ 3849566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 38533f72c5dSHong Zhang /* lambda_s^T F_PP w_2 */ 3869566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 387c9681e5dSHong Zhang } 388cd4cee2dSHong Zhang 389c9681e5dSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 3909566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 3919566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj])); 392cd4cee2dSHong Zhang if (quadJp) { 3939566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 3949566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 3959566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col)); 3969566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdp_col)); 3979566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 39833f72c5dSHong Zhang } 399c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 4009566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 4019566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj])); 40248a46eb9SPierre Jolivet if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj])); 40348a46eb9SPierre Jolivet if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj])); 404c9681e5dSHong Zhang } 405c9681e5dSHong Zhang } 406c9681e5dSHong Zhang } 407c9681e5dSHong Zhang 408c9681e5dSHong Zhang if (ts->vecs_sensi2) { 4099566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 4109566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 411c9681e5dSHong Zhang } 412c9681e5dSHong Zhang th->status = TS_STEP_COMPLETE; 4133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 414c9681e5dSHong Zhang } 415c9681e5dSHong Zhang 416d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_Theta(TS ts) 417d71ae5a4SJacob Faibussowitsch { 4182ca6e920SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 419cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 420b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp; 421b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp; 4222ca6e920SHong Zhang PetscInt nadj; 4237207e0fdSHong Zhang Mat J, Jpre, quadJ = NULL, quadJp = NULL; 4242ca6e920SHong Zhang KSP ksp; 425b5b4017aSHong Zhang PetscScalar *xarr; 4261a352fa8SHong Zhang PetscReal adjoint_time_step; 427da81f932SPierre 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) */ 4282ca6e920SHong Zhang 4292ca6e920SHong Zhang PetscFunctionBegin; 430c9681e5dSHong Zhang if (th->Theta == 1.) { 4319566063dSJacob Faibussowitsch PetscCall(TSAdjointStepBEuler_Private(ts)); 4323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 433c9681e5dSHong Zhang } 4342ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 4359566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(ts->snes, &ksp)); 4369566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 4377207e0fdSHong Zhang if (quadts) { 4389566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 4399566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 4407207e0fdSHong Zhang } 441b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 4421a352fa8SHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step); 4431a352fa8SHong Zhang adjoint_ptime = ts->ptime + ts->time_step; 4441a352fa8SHong Zhang adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */ 4452ca6e920SHong Zhang 44682ad9101SHong Zhang if (!th->endpoint) { 4475072d23cSHong Zhang /* recover th->X0 using vec_sol and the stage value th->X */ 4489566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol)); 44982ad9101SHong Zhang } 45082ad9101SHong Zhang 451b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 45233f72c5dSHong Zhang /* Cost function has an integral term */ 4537207e0fdSHong Zhang if (quadts) { 45405755b9cSHong Zhang if (th->endpoint) { 4559566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 45605755b9cSHong Zhang } else { 4579566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 45805755b9cSHong Zhang } 4597207e0fdSHong Zhang } 46033f72c5dSHong Zhang 461abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 4629566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); 4639566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step))); 464cd4cee2dSHong Zhang if (quadJ) { 4659566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 4669566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 4679566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col)); 4689566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdu_col)); 4699566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 47036eaed60SHong Zhang } 4712ca6e920SHong Zhang } 4723c54f92cSHong Zhang 473b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 4741a352fa8SHong Zhang th->shift = 1. / (th->Theta * adjoint_time_step); 4753c54f92cSHong Zhang if (th->endpoint) { 4769566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 4773c54f92cSHong Zhang } else { 4789566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); 4793c54f92cSHong Zhang } 4809566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 4812ca6e920SHong Zhang 482b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 483abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 4841d14d03bSHong Zhang KSPConvergedReason kspreason; 4859566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj])); 4869566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 4871d14d03bSHong Zhang if (kspreason < 0) { 4881d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 48963a3b9bcSJacob 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)); 4901d14d03bSHong Zhang } 4912ca6e920SHong Zhang } 4923c54f92cSHong Zhang 4931d14d03bSHong Zhang /* Second-order adjoint */ 494b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 4953c633725SBarry Smith PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta"); 496b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 4979566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 4989566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 499b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 5009566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 5019566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 5029566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 503b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 5049566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 505b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */ 5069566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 5079566063dSJacob Faibussowitsch PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift)); 5089566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj])); 50948a46eb9SPierre Jolivet if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj])); 510b5b4017aSHong Zhang } 511b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 512b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 5131d14d03bSHong Zhang KSPConvergedReason kspreason; 5149566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj])); 5159566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 5161d14d03bSHong Zhang if (kspreason < 0) { 5171d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 51863a3b9bcSJacob 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)); 5191d14d03bSHong Zhang } 520b5b4017aSHong Zhang } 521b5b4017aSHong Zhang } 522b5b4017aSHong Zhang 52336eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 5241d14d03bSHong Zhang if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 5251a352fa8SHong Zhang th->shift = 1. / ((th->Theta - 1.) * adjoint_time_step); 5261a352fa8SHong Zhang th->stage_time = adjoint_ptime; 5279566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre)); 5289566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 52933f72c5dSHong Zhang /* R_U at t_n */ 5301baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL)); 531abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 5329566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj])); 533cd4cee2dSHong Zhang if (quadJ) { 5349566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 5359566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 5369566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col)); 5379566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdu_col)); 5389566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 539b5b4017aSHong Zhang } 5409566063dSJacob Faibussowitsch PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift)); 541b5b4017aSHong Zhang } 5421d14d03bSHong Zhang 5431d14d03bSHong Zhang /* Second-order adjoint */ 5441d14d03bSHong Zhang if (ts->vecs_sensi2) { /* U_n */ 545b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 5469566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr)); 5479566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 548b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 5499566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 5509566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 5519566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr)); 552b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 5539566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 554b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 55533f72c5dSHong 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 */ 5569566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj])); 5579566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj])); 55848a46eb9SPierre Jolivet if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj])); 5599566063dSJacob Faibussowitsch PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift)); 560b5b4017aSHong Zhang } 561b5e0532dSHong Zhang } 5623fd52205SHong Zhang 563c586f2e8SHong Zhang th->stage_time = ts->ptime; /* recover the old value */ 564c586f2e8SHong Zhang 5653fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 566b5b4017aSHong Zhang /* U_{n+1} */ 56782ad9101SHong Zhang th->shift = 1.0 / (adjoint_time_step * th->Theta); 5689566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol)); 5699566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE)); 5701baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 571abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 5729566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 5739566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj])); 5743b711c3fSHong Zhang if (quadJp) { 5759566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 5769566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 5779566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col)); 5789566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdp_col)); 5799566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 5803b711c3fSHong Zhang } 5813fd52205SHong Zhang } 58233f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 58333f72c5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 5849566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 5859566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 586b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 5879566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 5889566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 5899566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 59033f72c5dSHong Zhang 591b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 5929566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 593b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 59433f72c5dSHong 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) */ 5959566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 5969566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj])); 59748a46eb9SPierre Jolivet if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj])); 59848a46eb9SPierre Jolivet if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj])); 599b5b4017aSHong Zhang } 600b5b4017aSHong Zhang } 601b5b4017aSHong Zhang 602b5b4017aSHong Zhang /* U_s */ 6039566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->Xdot)); 6049566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE)); 6051baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp)); 606abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 6079566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 6089566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj])); 6093b711c3fSHong Zhang if (quadJp) { 6109566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 6119566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 6129566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col)); 6139566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdp_col)); 6149566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 6153b711c3fSHong Zhang } 61633f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 61733f72c5dSHong Zhang /* Get w1 at t_n from TLM matrix */ 6189566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr)); 6199566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 620b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 6219566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 6229566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 6239566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr)); 624b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 6259566063dSJacob Faibussowitsch PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 626b5b4017aSHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 62733f72c5dSHong 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) */ 6289566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 6299566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj])); 63048a46eb9SPierre Jolivet if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj])); 63148a46eb9SPierre Jolivet if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj])); 63236eaed60SHong Zhang } 63336eaed60SHong Zhang } 6343fd52205SHong Zhang } 635b5e0532dSHong Zhang } 6363fd52205SHong Zhang } else { /* one-stage case */ 6373e05ceb1SHong Zhang th->shift = 0.0; 6389566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */ 6399566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 6401baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 641abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 6429566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj])); 6439566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj])); 644cd4cee2dSHong Zhang if (quadJ) { 6459566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 6469566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 6479566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col)); 6489566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdu_col)); 6499566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 65036eaed60SHong Zhang } 6512ca6e920SHong Zhang } 6523fd52205SHong Zhang if (ts->vecs_sensip) { 65382ad9101SHong Zhang th->shift = 1.0 / (adjoint_time_step * th->Theta); 6549566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); 6559566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 6561baa6e33SBarry Smith if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp)); 657abc2977eSBarry Smith for (nadj = 0; nadj < ts->numcost; nadj++) { 6589566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 6599566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj])); 660cd4cee2dSHong Zhang if (quadJp) { 6619566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 6629566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 6639566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col)); 6649566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_drdp_col)); 6659566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 66636eaed60SHong Zhang } 66736eaed60SHong Zhang } 6683fd52205SHong Zhang } 6692ca6e920SHong Zhang } 6702ca6e920SHong Zhang 6712ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6732ca6e920SHong Zhang } 6742ca6e920SHong Zhang 675d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X) 676d71ae5a4SJacob Faibussowitsch { 677cd652676SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 678be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 679cd652676SJed Brown 680cd652676SJed Brown PetscFunctionBegin; 6819566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X)); 682be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 6839566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, dt, th->Xdot, th->X)); 6843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 685cd652676SJed Brown } 686cd652676SJed Brown 687d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) 688d71ae5a4SJacob Faibussowitsch { 6891566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 6901566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 6911566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 6927453f775SEmil Constantinescu PetscReal wltea, wlter; 6931566a47fSLisandro Dalcin 6941566a47fSLisandro Dalcin PetscFunctionBegin; 6959371c9d4SSatish Balay if (!th->vec_sol_prev) { 6969371c9d4SSatish Balay *wlte = -1; 6973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6989371c9d4SSatish Balay } 6991566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 7009371c9d4SSatish Balay if (ts->steprestart) { 7019371c9d4SSatish Balay *wlte = -1; 7023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7039371c9d4SSatish Balay } 7041566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 705fecfb714SLisandro Dalcin { 706be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 707be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev / h; 7089371c9d4SSatish Balay PetscScalar scal[3]; 7099371c9d4SSatish Balay Vec vecs[3]; 71052f13ae7SStefano Zampini 711225f7d0bSStefano Zampini scal[0] = -1 / a; 712225f7d0bSStefano Zampini scal[1] = +1 / (a - 1); 713225f7d0bSStefano Zampini scal[2] = -1 / (a * (a - 1)); 7149371c9d4SSatish Balay vecs[0] = X; 7159371c9d4SSatish Balay vecs[1] = th->X0; 7169371c9d4SSatish Balay vecs[2] = th->vec_sol_prev; 7179566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Y)); 7189566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y, 3, scal, vecs)); 7199566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter)); 7201566a47fSLisandro Dalcin } 7211566a47fSLisandro Dalcin if (order) *order = 2; 7223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7231566a47fSLisandro Dalcin } 7241566a47fSLisandro Dalcin 725d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Theta(TS ts) 726d71ae5a4SJacob Faibussowitsch { 7271566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 7287207e0fdSHong Zhang TS quadts = ts->quadraturets; 7291566a47fSLisandro Dalcin 7301566a47fSLisandro Dalcin PetscFunctionBegin; 73148a46eb9SPierre Jolivet if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol)); 732715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 7331baa6e33SBarry Smith if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN)); 73448a46eb9SPierre Jolivet if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN)); 7353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 736715f1b00SHong Zhang } 737715f1b00SHong Zhang 738d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_Theta(TS ts) 739d71ae5a4SJacob Faibussowitsch { 740715f1b00SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 741cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 74213b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 743b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 7447207e0fdSHong Zhang PetscInt ntlm; 745715f1b00SHong Zhang KSP ksp; 7467207e0fdSHong Zhang Mat J, Jpre, quadJ = NULL, quadJp = NULL; 74713b1b0ffSHong Zhang PetscScalar *barr, *xarr; 7481a352fa8SHong Zhang PetscReal previous_shift; 749715f1b00SHong Zhang 750715f1b00SHong Zhang PetscFunctionBegin; 7511a352fa8SHong Zhang previous_shift = th->shift; 7529566063dSJacob Faibussowitsch PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN)); 75313b1b0ffSHong Zhang 75448a46eb9SPierre Jolivet if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN)); 7559566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(ts->snes, &ksp)); 7569566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 7577207e0fdSHong Zhang if (quadts) { 7589566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 7599566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 7607207e0fdSHong Zhang } 761715f1b00SHong Zhang 762715f1b00SHong Zhang /* Build RHS */ 763715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 7641a352fa8SHong Zhang th->shift = 1. / ((th->Theta - 1.) * th->time_step0); 7659566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 766fb842aefSJose E. Roman PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DETERMINE, &MatDeltaFwdSensip)); 7679566063dSJacob Faibussowitsch PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta)); 768715f1b00SHong Zhang 769022c081aSHong Zhang /* Add the f_p forcing terms */ 7700e11c21fSHong Zhang if (ts->Jacp) { 7719566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->Xdot)); 7729566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 7739566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN)); 77482ad9101SHong Zhang th->shift = previous_shift; 7759566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol)); 7769566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 7779566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN)); 7780e11c21fSHong Zhang } 779715f1b00SHong Zhang } else { /* 1-stage method */ 7803e05ceb1SHong Zhang th->shift = 0.0; 7819566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 782fb842aefSJose E. Roman PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DETERMINE, &MatDeltaFwdSensip)); 7839566063dSJacob Faibussowitsch PetscCall(MatScale(MatDeltaFwdSensip, -1.)); 78413b1b0ffSHong Zhang 785022c081aSHong Zhang /* Add the f_p forcing terms */ 7860e11c21fSHong Zhang if (ts->Jacp) { 78782ad9101SHong Zhang th->shift = previous_shift; 7889566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); 7899566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 7909566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN)); 791715f1b00SHong Zhang } 7920e11c21fSHong Zhang } 793715f1b00SHong Zhang 794715f1b00SHong Zhang /* Build LHS */ 7951a352fa8SHong Zhang th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */ 796715f1b00SHong Zhang if (th->endpoint) { 7979566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 798715f1b00SHong Zhang } else { 7999566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 800715f1b00SHong Zhang } 8019566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, J, Jpre)); 802715f1b00SHong Zhang 803715f1b00SHong Zhang /* 804715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 805c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 806715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 807715f1b00SHong Zhang */ 808715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 8097207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 8109566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL)); 8119566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp)); 812fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp)); 8139566063dSJacob Faibussowitsch PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 8149566063dSJacob Faibussowitsch PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 815715f1b00SHong Zhang } 816715f1b00SHong Zhang } 817715f1b00SHong Zhang 818715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 819022c081aSHong Zhang for (ntlm = 0; ntlm < th->num_tlm; ntlm++) { 820b5b4017aSHong Zhang KSPConvergedReason kspreason; 8219566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr)); 8229566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr)); 823715f1b00SHong Zhang if (th->endpoint) { 8249566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr)); 8259566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 8269566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col)); 8279566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 8289566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 829715f1b00SHong Zhang } else { 8309566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol)); 831715f1b00SHong Zhang } 8329566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 833b5b4017aSHong Zhang if (kspreason < 0) { 834b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 83563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm)); 836b5b4017aSHong Zhang } 8379566063dSJacob Faibussowitsch PetscCall(VecResetArray(VecDeltaFwdSensipCol)); 8389566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr)); 839715f1b00SHong Zhang } 840715f1b00SHong Zhang 84113b1b0ffSHong Zhang /* 84213b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 843c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 84413b1b0ffSHong Zhang */ 8457207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 84613b1b0ffSHong Zhang if (!th->endpoint) { 8479566063dSJacob Faibussowitsch PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */ 8489566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 8499566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp)); 850fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp)); 8519566063dSJacob Faibussowitsch PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 8529566063dSJacob Faibussowitsch PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 8539566063dSJacob Faibussowitsch PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); 85413b1b0ffSHong Zhang } else { 8559566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 8569566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 857fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp)); 8589566063dSJacob Faibussowitsch PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 8599566063dSJacob Faibussowitsch PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 86013b1b0ffSHong Zhang } 86113b1b0ffSHong Zhang } else { 86248a46eb9SPierre Jolivet if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); 863715f1b00SHong Zhang } 8643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8651566a47fSLisandro Dalcin } 8661566a47fSLisandro Dalcin 867d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[]) 868d71ae5a4SJacob Faibussowitsch { 8696affb6f8SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 8706affb6f8SHong Zhang 8716affb6f8SHong Zhang PetscFunctionBegin; 8727409eceaSHong Zhang if (ns) { 8737409eceaSHong Zhang if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */ 8747409eceaSHong Zhang else *ns = 2; /* endpoint form */ 8757409eceaSHong Zhang } 8765072d23cSHong Zhang if (stagesensip) { 877cc4f23bcSHong Zhang if (!th->endpoint && th->Theta != 1.0) { 8787409eceaSHong Zhang th->MatFwdStages[0] = th->MatDeltaFwdSensip; 879cc4f23bcSHong Zhang } else { 880cc4f23bcSHong Zhang th->MatFwdStages[0] = th->MatFwdSensip0; 881cc4f23bcSHong Zhang th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */ 882cc4f23bcSHong Zhang } 883cc4f23bcSHong Zhang *stagesensip = th->MatFwdStages; 8845072d23cSHong Zhang } 8853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8866affb6f8SHong Zhang } 8876affb6f8SHong Zhang 888316643e7SJed Brown /*------------------------------------------------------------*/ 889d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Theta(TS ts) 890d71ae5a4SJacob Faibussowitsch { 891316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 892316643e7SJed Brown 893316643e7SJed Brown PetscFunctionBegin; 8949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X)); 8959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Xdot)); 8969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X0)); 8979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->affine)); 8981566a47fSLisandro Dalcin 8999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_sol_prev)); 9009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work)); 9011566a47fSLisandro Dalcin 9029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->VecCostIntegral0)); 9033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 904277b19d0SLisandro Dalcin } 905277b19d0SLisandro Dalcin 906d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_Theta(TS ts) 907d71ae5a4SJacob Faibussowitsch { 908ece46509SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 909ece46509SHong Zhang 910ece46509SHong Zhang PetscFunctionBegin; 9119566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam)); 9129566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu)); 9139566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2)); 9149566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2)); 9159566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp)); 9169566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp)); 9173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 918ece46509SHong Zhang } 919ece46509SHong Zhang 920d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Theta(TS ts) 921d71ae5a4SJacob Faibussowitsch { 922277b19d0SLisandro Dalcin PetscFunctionBegin; 9239566063dSJacob Faibussowitsch PetscCall(TSReset_Theta(ts)); 924b3a6b972SJed Brown if (ts->dm) { 9259566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts)); 9269566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts)); 927b3a6b972SJed Brown } 9289566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 9299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL)); 9309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL)); 9319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL)); 9329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL)); 9333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 934316643e7SJed Brown } 935316643e7SJed Brown 936316643e7SJed Brown /* 937316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9382b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 9390fd10508SMatthew Knepley 9400fd10508SMatthew Knepley Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true, 9410fd10508SMatthew Knepley otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is 9420fd10508SMatthew Knepley U = (U_{n+1} + U0)/2 943316643e7SJed Brown */ 944d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts) 945d71ae5a4SJacob Faibussowitsch { 946316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 9477445fe48SJed Brown Vec X0, Xdot; 9487445fe48SJed Brown DM dm, dmsave; 9493e05ceb1SHong Zhang PetscReal shift = th->shift; 950316643e7SJed Brown 951316643e7SJed Brown PetscFunctionBegin; 9529566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 9535a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9549566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot)); 955494ce9fcSHong Zhang if (x != X0) { 9569566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); 957494ce9fcSHong Zhang } else { 9589566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Xdot)); 959494ce9fcSHong Zhang } 9607445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9617445fe48SJed Brown dmsave = ts->dm; 9627445fe48SJed Brown ts->dm = dm; 9639566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE)); 9647445fe48SJed Brown ts->dm = dmsave; 9659566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 9663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 967316643e7SJed Brown } 968316643e7SJed Brown 969d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts) 970d71ae5a4SJacob Faibussowitsch { 971316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 9727445fe48SJed Brown Vec Xdot; 9737445fe48SJed Brown DM dm, dmsave; 9743e05ceb1SHong Zhang PetscReal shift = th->shift; 975316643e7SJed Brown 976316643e7SJed Brown PetscFunctionBegin; 9779566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 978be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 9799566063dSJacob Faibussowitsch PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot)); 9807445fe48SJed Brown 9817445fe48SJed Brown dmsave = ts->dm; 9827445fe48SJed Brown ts->dm = dm; 9839566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE)); 9847445fe48SJed Brown ts->dm = dmsave; 9859566063dSJacob Faibussowitsch PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 9863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 987316643e7SJed Brown } 988316643e7SJed Brown 989d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_Theta(TS ts) 990d71ae5a4SJacob Faibussowitsch { 991715f1b00SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 9927207e0fdSHong Zhang TS quadts = ts->quadraturets; 993715f1b00SHong Zhang 994715f1b00SHong Zhang PetscFunctionBegin; 995022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 99613b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 9979566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip)); 9987207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 9999566063dSJacob Faibussowitsch PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp)); 10009566063dSJacob Faibussowitsch PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0)); 1001715f1b00SHong Zhang } 10026f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 10039566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0)); 100413b1b0ffSHong Zhang 10059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol)); 10063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1007715f1b00SHong Zhang } 1008715f1b00SHong Zhang 1009d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_Theta(TS ts) 1010d71ae5a4SJacob Faibussowitsch { 10117adebddeSHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 10127207e0fdSHong Zhang TS quadts = ts->quadraturets; 10137adebddeSHong Zhang 10147adebddeSHong Zhang PetscFunctionBegin; 10157207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 10169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&th->MatIntegralSensipTemp)); 10179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&th->MatIntegralSensip0)); 10187adebddeSHong Zhang } 10199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol)); 10209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&th->MatDeltaFwdSensip)); 10219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&th->MatFwdSensip0)); 10223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10237adebddeSHong Zhang } 10247adebddeSHong Zhang 1025d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Theta(TS ts) 1026d71ae5a4SJacob Faibussowitsch { 1027316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 1028cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 10292ffb9264SLisandro Dalcin PetscBool match; 1030316643e7SJed Brown 1031316643e7SJed Brown PetscFunctionBegin; 1032cd4cee2dSHong Zhang if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */ 10339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0)); 1034b1cb36f3SHong Zhang } 103548a46eb9SPierre Jolivet if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X)); 103648a46eb9SPierre Jolivet if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot)); 103748a46eb9SPierre Jolivet if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0)); 103848a46eb9SPierre Jolivet if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine)); 10393b1890cdSShri Abhyankar 10401566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 104120d49056SMatthew G. Knepley th->shift = 1 / (th->Theta * ts->time_step); 10421566a47fSLisandro Dalcin 10439566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &ts->dm)); 10449566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts)); 10459566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts)); 10461566a47fSLisandro Dalcin 10479566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 10489566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 10499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match)); 10502ffb9264SLisandro Dalcin if (!match) { 1051ecc87898SStefano Zampini if (!th->vec_sol_prev) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev)); 1052ecc87898SStefano Zampini if (!th->vec_lte_work) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work)); 10533b1890cdSShri Abhyankar } 10549566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &ts->snes)); 1055cc4f23bcSHong Zhang 1056cc4f23bcSHong Zhang ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE; 10573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1058b4dd3bf9SBarry Smith } 10590c7376e5SHong Zhang 1060b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1061b4dd3bf9SBarry Smith 1062d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1063d71ae5a4SJacob Faibussowitsch { 1064b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta *)ts->data; 1065b4dd3bf9SBarry Smith 1066b4dd3bf9SBarry Smith PetscFunctionBegin; 10679566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam)); 10689566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp)); 106948a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu)); 1070b5b4017aSHong Zhang if (ts->vecs_sensi2) { 10719566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2)); 10729566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp)); 107367633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 107467633408SHong Zhang if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu; 107567633408SHong Zhang if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup; 1076b5b4017aSHong Zhang } 1077c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 10789566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2)); 107967633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 108067633408SHong Zhang if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu; 108167633408SHong Zhang if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp; 1082b5b4017aSHong Zhang } 10833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1084316643e7SJed Brown } 1085316643e7SJed Brown 1086ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems PetscOptionsObject) 1087d71ae5a4SJacob Faibussowitsch { 1088316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 1089316643e7SJed Brown 1090316643e7SJed Brown PetscFunctionBegin; 1091d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options"); 1092316643e7SJed Brown { 10939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL)); 10949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL)); 10959566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL)); 1096316643e7SJed Brown } 1097d0609cedSBarry Smith PetscOptionsHeadEnd(); 10983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1099316643e7SJed Brown } 1100316643e7SJed Brown 1101d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer) 1102d71ae5a4SJacob Faibussowitsch { 1103316643e7SJed Brown TS_Theta *th = (TS_Theta *)ts->data; 11049f196a02SMartin Diehl PetscBool isascii; 1105316643e7SJed Brown 1106316643e7SJed Brown PetscFunctionBegin; 11079f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 11089f196a02SMartin Diehl if (isascii) { 11099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Theta=%g\n", (double)th->Theta)); 11109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Extrapolation=%s\n", th->extrapolate ? "yes" : "no")); 1111316643e7SJed Brown } 11123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1113316643e7SJed Brown } 1114316643e7SJed Brown 1115d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta) 1116d71ae5a4SJacob Faibussowitsch { 11170de4c49aSJed Brown TS_Theta *th = (TS_Theta *)ts->data; 11180de4c49aSJed Brown 11190de4c49aSJed Brown PetscFunctionBegin; 11200de4c49aSJed Brown *theta = th->Theta; 11213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11220de4c49aSJed Brown } 11230de4c49aSJed Brown 1124d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta) 1125d71ae5a4SJacob Faibussowitsch { 11260de4c49aSJed Brown TS_Theta *th = (TS_Theta *)ts->data; 11270de4c49aSJed Brown 11280de4c49aSJed Brown PetscFunctionBegin; 1129cad9d221SBarry Smith PetscCheck(theta > 0 && theta <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Theta %g not in range (0,1]", (double)theta); 11300de4c49aSJed Brown th->Theta = theta; 11311566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 11323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11330de4c49aSJed Brown } 1134eb284becSJed Brown 1135d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint) 1136d71ae5a4SJacob Faibussowitsch { 113726f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 113826f2ff8fSLisandro Dalcin 113926f2ff8fSLisandro Dalcin PetscFunctionBegin; 114026f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 11413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114226f2ff8fSLisandro Dalcin } 114326f2ff8fSLisandro Dalcin 1144d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg) 1145d71ae5a4SJacob Faibussowitsch { 1146eb284becSJed Brown TS_Theta *th = (TS_Theta *)ts->data; 1147eb284becSJed Brown 1148eb284becSJed Brown PetscFunctionBegin; 1149eb284becSJed Brown th->endpoint = flg; 11503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1151eb284becSJed Brown } 11520de4c49aSJed Brown 1153f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1154d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) 1155d71ae5a4SJacob Faibussowitsch { 1156f9c1d6abSBarry Smith PetscComplex z = xr + xi * PETSC_i, f; 1157f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta *)ts->data; 1158f9c1d6abSBarry Smith 1159f9c1d6abSBarry Smith PetscFunctionBegin; 1160520136c7SJose E. Roman f = (1.0 + (1.0 - th->Theta) * z) / (1.0 - th->Theta * z); 1161f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 1162f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 11633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1164f9c1d6abSBarry Smith } 1165f9c1d6abSBarry Smith #endif 1166f9c1d6abSBarry Smith 1167d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[]) 1168d71ae5a4SJacob Faibussowitsch { 116942682096SHong Zhang TS_Theta *th = (TS_Theta *)ts->data; 117042682096SHong Zhang 117142682096SHong Zhang PetscFunctionBegin; 11727409eceaSHong Zhang if (ns) { 11737409eceaSHong Zhang if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */ 11747409eceaSHong Zhang else *ns = 2; /* endpoint form */ 11757409eceaSHong Zhang } 11765072d23cSHong Zhang if (Y) { 1177cc4f23bcSHong Zhang if (!th->endpoint && th->Theta != 1.0) { 11787409eceaSHong Zhang th->Stages[0] = th->X; 1179cc4f23bcSHong Zhang } else { 1180cc4f23bcSHong Zhang th->Stages[0] = th->X0; 1181cc4f23bcSHong Zhang th->Stages[1] = ts->vec_sol; /* stiffly accurate */ 1182cc4f23bcSHong Zhang } 1183cc4f23bcSHong Zhang *Y = th->Stages; 11845072d23cSHong Zhang } 11853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118642682096SHong Zhang } 1187f9c1d6abSBarry Smith 1188316643e7SJed Brown /* ------------------------------------------------------------ */ 1189316643e7SJed Brown /*MC 119096f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1191316643e7SJed Brown 1192316643e7SJed Brown Level: beginner 1193316643e7SJed Brown 1194bcf0153eSBarry Smith Options Database Keys: 1195c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1196c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 119703842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 11984eb428fdSBarry Smith 1199eb284becSJed Brown Notes: 120037fdd005SBarry Smith .vb 120137fdd005SBarry Smith -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 120237fdd005SBarry Smith -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 120337fdd005SBarry Smith -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 120437fdd005SBarry Smith .ve 12054eb428fdSBarry Smith 12067409eceaSHong 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. 1207eb284becSJed Brown 12087409eceaSHong Zhang The midpoint variant is cast as a 1-stage implicit Runge-Kutta method. 1209eb284becSJed Brown 1210eb284becSJed Brown .vb 1211eb284becSJed Brown Theta | Theta 1212eb284becSJed Brown ------------- 1213eb284becSJed Brown | 1 1214eb284becSJed Brown .ve 1215eb284becSJed Brown 1216eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1217eb284becSJed Brown 1218eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1219eb284becSJed Brown 1220eb284becSJed Brown .vb 1221eb284becSJed Brown 0 | 0 0 1222eb284becSJed Brown 1 | 1-Theta Theta 1223eb284becSJed Brown ------------------- 1224eb284becSJed Brown | 1-Theta Theta 1225eb284becSJed Brown .ve 1226eb284becSJed Brown 1227eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1228eb284becSJed Brown 1229eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1230b44f4de4SBarry Smith .vb 1231b44f4de4SBarry Smith Y_i = X + h sum_j a_ij Y'_j 1232b44f4de4SBarry Smith .ve 12334eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1234eb284becSJed Brown 12351cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()` 1236316643e7SJed Brown M*/ 1237d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1238d71ae5a4SJacob Faibussowitsch { 1239316643e7SJed Brown TS_Theta *th; 1240316643e7SJed Brown 1241316643e7SJed Brown PetscFunctionBegin; 1242277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1243ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1244316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1245316643e7SJed Brown ts->ops->view = TSView_Theta; 1246316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 124742f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1248ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1249316643e7SJed Brown ts->ops->step = TSStep_Theta; 1250cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12511566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 125224655328SShri ts->ops->rollback = TSRollBack_Theta; 125326b5f05dSStefano Zampini ts->ops->resizeregister = TSResizeRegister_Theta; 1254316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12550f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12560f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1257f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1258f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1259f9c1d6abSBarry Smith #endif 126042682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 126142f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1262b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1263b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12642ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12656affb6f8SHong Zhang 1266715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12677adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1268715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12696affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1270316643e7SJed Brown 1271efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1272efd4aadfSBarry Smith 12734dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&th)); 1274316643e7SJed Brown ts->data = (void *)th; 1275316643e7SJed Brown 1276715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1277715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1278715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1279b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1280715f1b00SHong Zhang 12816f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1282316643e7SJed Brown th->Theta = 0.5; 12831566a47fSLisandro Dalcin th->order = 2; 12849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta)); 12859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta)); 12869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta)); 12879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta)); 12883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1289316643e7SJed Brown } 12900de4c49aSJed Brown 12910de4c49aSJed Brown /*@ 1292bcf0153eSBarry Smith TSThetaGetTheta - Get the abscissa of the stage in (0,1] for `TSTHETA` 12930de4c49aSJed Brown 12940de4c49aSJed Brown Not Collective 12950de4c49aSJed Brown 12960de4c49aSJed Brown Input Parameter: 12970de4c49aSJed Brown . ts - timestepping context 12980de4c49aSJed Brown 12990de4c49aSJed Brown Output Parameter: 13000de4c49aSJed Brown . theta - stage abscissa 13010de4c49aSJed Brown 1302bcf0153eSBarry Smith Level: advanced 1303bcf0153eSBarry Smith 13040de4c49aSJed Brown Note: 1305bcf0153eSBarry Smith Use of this function is normally only required to hack `TSTHETA` to use a modified integration scheme. 13060de4c49aSJed Brown 13071cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetTheta()`, `TSTHETA` 13080de4c49aSJed Brown @*/ 1309d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta) 1310d71ae5a4SJacob Faibussowitsch { 13110de4c49aSJed Brown PetscFunctionBegin; 1312afb20b64SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 13134f572ea9SToby Isaac PetscAssertPointer(theta, 2); 1314cac4c232SBarry Smith PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta)); 13153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13160de4c49aSJed Brown } 13170de4c49aSJed Brown 13180de4c49aSJed Brown /*@ 1319bcf0153eSBarry Smith TSThetaSetTheta - Set the abscissa of the stage in (0,1] for `TSTHETA` 13200de4c49aSJed Brown 13210de4c49aSJed Brown Not Collective 13220de4c49aSJed Brown 1323d8d19677SJose E. Roman Input Parameters: 13240de4c49aSJed Brown + ts - timestepping context 13250de4c49aSJed Brown - theta - stage abscissa 13260de4c49aSJed Brown 1327bcf0153eSBarry Smith Options Database Key: 132867b8a455SSatish Balay . -ts_theta_theta <theta> - set theta 13290de4c49aSJed Brown 1330bcf0153eSBarry Smith Level: intermediate 13310de4c49aSJed Brown 13321cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaGetTheta()`, `TSTHETA`, `TSCN` 13330de4c49aSJed Brown @*/ 1334d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta) 1335d71ae5a4SJacob Faibussowitsch { 13360de4c49aSJed Brown PetscFunctionBegin; 1337afb20b64SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1338cac4c232SBarry Smith PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta)); 13393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13400de4c49aSJed Brown } 1341f33bbcb6SJed Brown 134226f2ff8fSLisandro Dalcin /*@ 1343bcf0153eSBarry Smith TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA` 134426f2ff8fSLisandro Dalcin 134526f2ff8fSLisandro Dalcin Not Collective 134626f2ff8fSLisandro Dalcin 134726f2ff8fSLisandro Dalcin Input Parameter: 134826f2ff8fSLisandro Dalcin . ts - timestepping context 134926f2ff8fSLisandro Dalcin 135026f2ff8fSLisandro Dalcin Output Parameter: 1351bcf0153eSBarry Smith . endpoint - `PETSC_TRUE` when using the endpoint variant 135226f2ff8fSLisandro Dalcin 1353bcf0153eSBarry Smith Level: advanced 135426f2ff8fSLisandro Dalcin 13551cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN` 135626f2ff8fSLisandro Dalcin @*/ 1357d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint) 1358d71ae5a4SJacob Faibussowitsch { 135926f2ff8fSLisandro Dalcin PetscFunctionBegin; 136026f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 13614f572ea9SToby Isaac PetscAssertPointer(endpoint, 2); 1362cac4c232SBarry Smith PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint)); 13633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136426f2ff8fSLisandro Dalcin } 136526f2ff8fSLisandro Dalcin 1366eb284becSJed Brown /*@ 1367bcf0153eSBarry Smith TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA` 1368eb284becSJed Brown 1369eb284becSJed Brown Not Collective 1370eb284becSJed Brown 1371d8d19677SJose E. Roman Input Parameters: 1372eb284becSJed Brown + ts - timestepping context 1373bcf0153eSBarry Smith - flg - `PETSC_TRUE` to use the endpoint variant 1374eb284becSJed Brown 1375bcf0153eSBarry Smith Options Database Key: 137667b8a455SSatish Balay . -ts_theta_endpoint <flg> - use the endpoint variant 1377eb284becSJed Brown 1378bcf0153eSBarry Smith Level: intermediate 1379eb284becSJed Brown 13801cc06b55SBarry Smith .seealso: [](ch_ts), `TSTHETA`, `TSCN` 1381eb284becSJed Brown @*/ 1382d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg) 1383d71ae5a4SJacob Faibussowitsch { 1384eb284becSJed Brown PetscFunctionBegin; 1385eb284becSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1386cac4c232SBarry Smith PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg)); 13873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1388eb284becSJed Brown } 1389eb284becSJed Brown 1390f33bbcb6SJed Brown /* 1391f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1392f33bbcb6SJed Brown * The creation functions for these specializations are below. 1393f33bbcb6SJed Brown */ 1394f33bbcb6SJed Brown 1395d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BEuler(TS ts) 1396d71ae5a4SJacob Faibussowitsch { 13971566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 13981566a47fSLisandro Dalcin 13991566a47fSLisandro Dalcin PetscFunctionBegin; 140008401ef6SPierre 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"); 140128b400f6SJacob 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"); 14029566063dSJacob Faibussowitsch PetscCall(TSSetUp_Theta(ts)); 14033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14041566a47fSLisandro Dalcin } 14051566a47fSLisandro Dalcin 1406d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer) 1407d71ae5a4SJacob Faibussowitsch { 1408f33bbcb6SJed Brown PetscFunctionBegin; 14093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1410f33bbcb6SJed Brown } 1411f33bbcb6SJed Brown 1412f33bbcb6SJed Brown /*MC 1413f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1414f33bbcb6SJed Brown 1415f33bbcb6SJed Brown Level: beginner 1416f33bbcb6SJed Brown 1417bcf0153eSBarry Smith Note: 141837fdd005SBarry Smith `TSBEULER` is equivalent to `TSTHETA` with Theta=1.0 or `-ts_type theta -ts_theta_theta 1.0` 14194eb428fdSBarry Smith 14201cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA` 1421f33bbcb6SJed Brown M*/ 1422d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1423d71ae5a4SJacob Faibussowitsch { 1424f33bbcb6SJed Brown PetscFunctionBegin; 14259566063dSJacob Faibussowitsch PetscCall(TSCreate_Theta(ts)); 14269566063dSJacob Faibussowitsch PetscCall(TSThetaSetTheta(ts, 1.0)); 14279566063dSJacob Faibussowitsch PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE)); 14280c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1429f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 14303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1431f33bbcb6SJed Brown } 1432f33bbcb6SJed Brown 1433d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_CN(TS ts) 1434d71ae5a4SJacob Faibussowitsch { 14351566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta *)ts->data; 14361566a47fSLisandro Dalcin 14371566a47fSLisandro Dalcin PetscFunctionBegin; 143808401ef6SPierre 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"); 14393c633725SBarry 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"); 14409566063dSJacob Faibussowitsch PetscCall(TSSetUp_Theta(ts)); 14413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14421566a47fSLisandro Dalcin } 14431566a47fSLisandro Dalcin 1444d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer) 1445d71ae5a4SJacob Faibussowitsch { 1446f33bbcb6SJed Brown PetscFunctionBegin; 14473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1448f33bbcb6SJed Brown } 1449f33bbcb6SJed Brown 1450f33bbcb6SJed Brown /*MC 1451f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1452f33bbcb6SJed Brown 1453f33bbcb6SJed Brown Level: beginner 1454f33bbcb6SJed Brown 1455f33bbcb6SJed Brown Notes: 1456bcf0153eSBarry Smith `TSCN` is equivalent to `TSTHETA` with Theta=0.5 and the "endpoint" option set. I.e. 1457bcf0153eSBarry Smith .vb 1458bcf0153eSBarry Smith -ts_type theta 1459bcf0153eSBarry Smith -ts_theta_theta 0.5 1460bcf0153eSBarry Smith -ts_theta_endpoint 1461bcf0153eSBarry Smith .ve 14627cf5af47SJed Brown 14631cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`, `TSType`, 1464f33bbcb6SJed Brown M*/ 1465d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1466d71ae5a4SJacob Faibussowitsch { 1467f33bbcb6SJed Brown PetscFunctionBegin; 14689566063dSJacob Faibussowitsch PetscCall(TSCreate_Theta(ts)); 14699566063dSJacob Faibussowitsch PetscCall(TSThetaSetTheta(ts, 0.5)); 14709566063dSJacob Faibussowitsch PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE)); 14710c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1472f33bbcb6SJed Brown ts->ops->view = TSView_CN; 14733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1474f33bbcb6SJed Brown } 1475