1 /* 2 Code for timestepping with implicit Theta method 3 */ 4 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5 #include <petscsnes.h> 6 #include <petscdm.h> 7 #include <petscmat.h> 8 9 typedef struct { 10 /* context for time stepping */ 11 PetscReal stage_time; 12 Vec Stages[2]; /* Storage for stage solutions */ 13 Vec X0, X, Xdot; /* Storage for u^n, u^n + dt a_{11} k_1, and time derivative u^{n+1}_t */ 14 Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 15 PetscReal Theta; 16 PetscReal shift; /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */ 17 PetscInt order; 18 PetscBool endpoint; 19 PetscBool extrapolate; 20 TSStepStatus status; 21 Vec VecCostIntegral0; /* Backup for roll-backs due to events, used by cost integral */ 22 PetscReal ptime0; /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */ 23 PetscReal time_step0; /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/ 24 25 /* context for sensitivity analysis */ 26 PetscInt num_tlm; /* Total number of tangent linear equations */ 27 Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 28 Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 29 Vec *VecsSensiTemp; /* Vector to be multiplied with Jacobian transpose */ 30 Mat MatFwdStages[2]; /* TLM Stages */ 31 Mat MatDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 32 Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */ 33 Mat MatFwdSensip0; /* backup for roll-backs due to events */ 34 Mat MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 35 Mat MatIntegralSensip0; /* backup for roll-backs due to events */ 36 Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */ 37 Vec *VecsDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */ 38 Vec *VecsSensi2Temp; /* Working vectors that holds the residual for the second-order adjoint */ 39 Vec *VecsAffine; /* Working vectors to store residuals */ 40 /* context for error estimation */ 41 Vec vec_sol_prev; 42 Vec vec_lte_work; 43 } TS_Theta; 44 45 static PetscErrorCode TSThetaGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 46 { 47 TS_Theta *th = (TS_Theta *)ts->data; 48 49 PetscFunctionBegin; 50 if (X0) { 51 if (dm && dm != ts->dm) { 52 PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_X0", X0)); 53 } else *X0 = ts->vec_sol; 54 } 55 if (Xdot) { 56 if (dm && dm != ts->dm) { 57 PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_Xdot", Xdot)); 58 } else *Xdot = th->Xdot; 59 } 60 PetscFunctionReturn(PETSC_SUCCESS); 61 } 62 63 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 64 { 65 PetscFunctionBegin; 66 if (X0) { 67 if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_X0", X0)); 68 } 69 if (Xdot) { 70 if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_Xdot", Xdot)); 71 } 72 PetscFunctionReturn(PETSC_SUCCESS); 73 } 74 75 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine, DM coarse, void *ctx) 76 { 77 PetscFunctionBegin; 78 PetscFunctionReturn(PETSC_SUCCESS); 79 } 80 81 static PetscErrorCode DMRestrictHook_TSTheta(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 82 { 83 TS ts = (TS)ctx; 84 Vec X0, Xdot, X0_c, Xdot_c; 85 86 PetscFunctionBegin; 87 PetscCall(TSThetaGetX0AndXdot(ts, fine, &X0, &Xdot)); 88 PetscCall(TSThetaGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 89 PetscCall(MatRestrict(restrct, X0, X0_c)); 90 PetscCall(MatRestrict(restrct, Xdot, Xdot_c)); 91 PetscCall(VecPointwiseMult(X0_c, rscale, X0_c)); 92 PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c)); 93 PetscCall(TSThetaRestoreX0AndXdot(ts, fine, &X0, &Xdot)); 94 PetscCall(TSThetaRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 95 PetscFunctionReturn(PETSC_SUCCESS); 96 } 97 98 static PetscErrorCode DMSubDomainHook_TSTheta(DM dm, DM subdm, void *ctx) 99 { 100 PetscFunctionBegin; 101 PetscFunctionReturn(PETSC_SUCCESS); 102 } 103 104 static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 105 { 106 TS ts = (TS)ctx; 107 Vec X0, Xdot, X0_sub, Xdot_sub; 108 109 PetscFunctionBegin; 110 PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot)); 111 PetscCall(TSThetaGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 112 113 PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 114 PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 115 116 PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 117 PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 118 119 PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 120 PetscCall(TSThetaRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 121 PetscFunctionReturn(PETSC_SUCCESS); 122 } 123 124 static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 125 { 126 TS_Theta *th = (TS_Theta *)ts->data; 127 TS quadts = ts->quadraturets; 128 129 PetscFunctionBegin; 130 if (th->endpoint) { 131 /* Evolve ts->vec_costintegral to compute integrals */ 132 if (th->Theta != 1.0) { 133 PetscCall(TSComputeRHSFunction(quadts, th->ptime0, th->X0, ts->vec_costintegrand)); 134 PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * (1.0 - th->Theta), ts->vec_costintegrand)); 135 } 136 PetscCall(TSComputeRHSFunction(quadts, ts->ptime, ts->vec_sol, ts->vec_costintegrand)); 137 PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * th->Theta, ts->vec_costintegrand)); 138 } else { 139 PetscCall(TSComputeRHSFunction(quadts, th->stage_time, th->X, ts->vec_costintegrand)); 140 PetscCall(VecAXPY(quadts->vec_sol, th->time_step0, ts->vec_costintegrand)); 141 } 142 PetscFunctionReturn(PETSC_SUCCESS); 143 } 144 145 static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 146 { 147 TS_Theta *th = (TS_Theta *)ts->data; 148 TS quadts = ts->quadraturets; 149 150 PetscFunctionBegin; 151 /* backup cost integral */ 152 PetscCall(VecCopy(quadts->vec_sol, th->VecCostIntegral0)); 153 PetscCall(TSThetaEvaluateCostIntegral(ts)); 154 PetscFunctionReturn(PETSC_SUCCESS); 155 } 156 157 static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 158 { 159 TS_Theta *th = (TS_Theta *)ts->data; 160 161 PetscFunctionBegin; 162 /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */ 163 th->ptime0 = ts->ptime + ts->time_step; 164 th->time_step0 = -ts->time_step; 165 PetscCall(TSThetaEvaluateCostIntegral(ts)); 166 PetscFunctionReturn(PETSC_SUCCESS); 167 } 168 169 static PetscErrorCode TSTheta_SNESSolve(TS ts, Vec b, Vec x) 170 { 171 PetscInt nits, lits; 172 173 PetscFunctionBegin; 174 PetscCall(SNESSolve(ts->snes, b, x)); 175 PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 176 PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 177 ts->snes_its += nits; 178 ts->ksp_its += lits; 179 PetscFunctionReturn(PETSC_SUCCESS); 180 } 181 182 /* We need to transfer X0 which will be copied into sol_prev */ 183 static PetscErrorCode TSResizeRegister_Theta(TS ts, PetscBool reg) 184 { 185 TS_Theta *th = (TS_Theta *)ts->data; 186 const char name[] = "ts:theta:X0"; 187 188 PetscFunctionBegin; 189 if (reg && th->vec_sol_prev) { 190 PetscCall(TSResizeRegisterVec(ts, name, th->X0)); 191 } else if (!reg) { 192 PetscCall(TSResizeRetrieveVec(ts, name, &th->X0)); 193 PetscCall(PetscObjectReference((PetscObject)th->X0)); 194 } 195 PetscFunctionReturn(PETSC_SUCCESS); 196 } 197 198 static PetscErrorCode TSStep_Theta(TS ts) 199 { 200 TS_Theta *th = (TS_Theta *)ts->data; 201 PetscInt rejections = 0; 202 PetscBool stageok, accept = PETSC_TRUE; 203 PetscReal next_time_step = ts->time_step; 204 205 PetscFunctionBegin; 206 if (!ts->steprollback) { 207 if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev)); 208 PetscCall(VecCopy(ts->vec_sol, th->X0)); 209 } 210 211 th->status = TS_STEP_INCOMPLETE; 212 while (!ts->reason && th->status != TS_STEP_COMPLETE) { 213 th->shift = 1 / (th->Theta * ts->time_step); 214 th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta) * ts->time_step; 215 PetscCall(VecCopy(th->X0, th->X)); 216 if (th->extrapolate && !ts->steprestart) PetscCall(VecAXPY(th->X, 1 / th->shift, th->Xdot)); 217 if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 218 if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol, &th->affine)); 219 PetscCall(VecZeroEntries(th->Xdot)); 220 PetscCall(TSComputeIFunction(ts, ts->ptime, th->X0, th->Xdot, th->affine, PETSC_FALSE)); 221 PetscCall(VecScale(th->affine, (th->Theta - 1) / th->Theta)); 222 } 223 PetscCall(TSPreStage(ts, th->stage_time)); 224 PetscCall(TSTheta_SNESSolve(ts, th->endpoint ? th->affine : NULL, th->X)); 225 PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X)); 226 PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok)); 227 if (!stageok) goto reject_step; 228 229 th->status = TS_STEP_PENDING; 230 if (th->endpoint) { 231 PetscCall(VecCopy(th->X, ts->vec_sol)); 232 } else { 233 PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); 234 PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot)); 235 } 236 PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 237 th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 238 if (!accept) { 239 PetscCall(VecCopy(th->X0, ts->vec_sol)); 240 ts->time_step = next_time_step; 241 goto reject_step; 242 } 243 244 if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 245 th->ptime0 = ts->ptime; 246 th->time_step0 = ts->time_step; 247 } 248 ts->ptime += ts->time_step; 249 ts->time_step = next_time_step; 250 break; 251 252 reject_step: 253 ts->reject++; 254 accept = PETSC_FALSE; 255 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 256 ts->reason = TS_DIVERGED_STEP_REJECTED; 257 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 258 } 259 } 260 PetscFunctionReturn(PETSC_SUCCESS); 261 } 262 263 static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 264 { 265 TS_Theta *th = (TS_Theta *)ts->data; 266 TS quadts = ts->quadraturets; 267 Vec *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp; 268 Vec *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp; 269 PetscInt nadj; 270 Mat J, Jpre, quadJ = NULL, quadJp = NULL; 271 KSP ksp; 272 PetscScalar *xarr; 273 TSEquationType eqtype; 274 PetscBool isexplicitode = PETSC_FALSE; 275 PetscReal adjoint_time_step; 276 277 PetscFunctionBegin; 278 PetscCall(TSGetEquationType(ts, &eqtype)); 279 if (eqtype == TS_EQ_ODE_EXPLICIT) { 280 isexplicitode = PETSC_TRUE; 281 VecsDeltaLam = ts->vecs_sensi; 282 VecsDeltaLam2 = ts->vecs_sensi2; 283 } 284 th->status = TS_STEP_INCOMPLETE; 285 PetscCall(SNESGetKSP(ts->snes, &ksp)); 286 PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 287 if (quadts) { 288 PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 289 PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 290 } 291 292 th->stage_time = ts->ptime; 293 adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */ 294 295 /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */ 296 if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 297 298 for (nadj = 0; nadj < ts->numcost; nadj++) { 299 PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); 300 PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */ 301 if (quadJ) { 302 PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 303 PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 304 PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col)); 305 PetscCall(VecResetArray(ts->vec_drdu_col)); 306 PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 307 } 308 } 309 310 /* Build LHS for first-order adjoint */ 311 th->shift = 1. / adjoint_time_step; 312 PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 313 PetscCall(KSPSetOperators(ksp, J, Jpre)); 314 315 /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 316 for (nadj = 0; nadj < ts->numcost; nadj++) { 317 KSPConvergedReason kspreason; 318 PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj])); 319 PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 320 if (kspreason < 0) { 321 ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 322 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj)); 323 } 324 } 325 326 if (ts->vecs_sensi2) { /* U_{n+1} */ 327 /* Get w1 at t_{n+1} from TLM matrix */ 328 PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 329 PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 330 /* lambda_s^T F_UU w_1 */ 331 PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 332 /* lambda_s^T F_UP w_2 */ 333 PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 334 for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */ 335 PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 336 PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step)); 337 PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj])); 338 if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj])); 339 } 340 /* Solve stage equation LHS X = RHS for second-order adjoint */ 341 for (nadj = 0; nadj < ts->numcost; nadj++) { 342 KSPConvergedReason kspreason; 343 PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj])); 344 PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 345 if (kspreason < 0) { 346 ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 347 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj)); 348 } 349 } 350 } 351 352 /* Update sensitivities, and evaluate integrals if there is any */ 353 if (!isexplicitode) { 354 th->shift = 0.0; 355 PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 356 PetscCall(KSPSetOperators(ksp, J, Jpre)); 357 for (nadj = 0; nadj < ts->numcost; nadj++) { 358 /* Add f_U \lambda_s to the original RHS */ 359 PetscCall(VecScale(VecsSensiTemp[nadj], -1.)); 360 PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj])); 361 PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step)); 362 PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj])); 363 if (ts->vecs_sensi2) { 364 PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj])); 365 PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step)); 366 PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj])); 367 } 368 } 369 } 370 if (ts->vecs_sensip) { 371 PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol)); 372 PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */ 373 if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 374 if (ts->vecs_sensi2p) { 375 /* lambda_s^T F_PU w_1 */ 376 PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 377 /* lambda_s^T F_PP w_2 */ 378 PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 379 } 380 381 for (nadj = 0; nadj < ts->numcost; nadj++) { 382 PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 383 PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj])); 384 if (quadJp) { 385 PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 386 PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 387 PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col)); 388 PetscCall(VecResetArray(ts->vec_drdp_col)); 389 PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 390 } 391 if (ts->vecs_sensi2p) { 392 PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 393 PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj])); 394 if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj])); 395 if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj])); 396 } 397 } 398 } 399 400 if (ts->vecs_sensi2) { 401 PetscCall(VecResetArray(ts->vec_sensip_col)); 402 PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 403 } 404 th->status = TS_STEP_COMPLETE; 405 PetscFunctionReturn(PETSC_SUCCESS); 406 } 407 408 static PetscErrorCode TSAdjointStep_Theta(TS ts) 409 { 410 TS_Theta *th = (TS_Theta *)ts->data; 411 TS quadts = ts->quadraturets; 412 Vec *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp; 413 Vec *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp; 414 PetscInt nadj; 415 Mat J, Jpre, quadJ = NULL, quadJp = NULL; 416 KSP ksp; 417 PetscScalar *xarr; 418 PetscReal adjoint_time_step; 419 PetscReal adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, usually ts->ptime is larger than adjoint_ptime) */ 420 421 PetscFunctionBegin; 422 if (th->Theta == 1.) { 423 PetscCall(TSAdjointStepBEuler_Private(ts)); 424 PetscFunctionReturn(PETSC_SUCCESS); 425 } 426 th->status = TS_STEP_INCOMPLETE; 427 PetscCall(SNESGetKSP(ts->snes, &ksp)); 428 PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 429 if (quadts) { 430 PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 431 PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 432 } 433 /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 434 th->stage_time = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step); 435 adjoint_ptime = ts->ptime + ts->time_step; 436 adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */ 437 438 if (!th->endpoint) { 439 /* recover th->X0 using vec_sol and the stage value th->X */ 440 PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol)); 441 } 442 443 /* Build RHS for first-order adjoint */ 444 /* Cost function has an integral term */ 445 if (quadts) { 446 if (th->endpoint) { 447 PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 448 } else { 449 PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 450 } 451 } 452 453 for (nadj = 0; nadj < ts->numcost; nadj++) { 454 PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); 455 PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step))); 456 if (quadJ) { 457 PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 458 PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 459 PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col)); 460 PetscCall(VecResetArray(ts->vec_drdu_col)); 461 PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 462 } 463 } 464 465 /* Build LHS for first-order adjoint */ 466 th->shift = 1. / (th->Theta * adjoint_time_step); 467 if (th->endpoint) { 468 PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre)); 469 } else { 470 PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); 471 } 472 PetscCall(KSPSetOperators(ksp, J, Jpre)); 473 474 /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 475 for (nadj = 0; nadj < ts->numcost; nadj++) { 476 KSPConvergedReason kspreason; 477 PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj])); 478 PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 479 if (kspreason < 0) { 480 ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 481 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj)); 482 } 483 } 484 485 /* Second-order adjoint */ 486 if (ts->vecs_sensi2) { /* U_{n+1} */ 487 PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta"); 488 /* Get w1 at t_{n+1} from TLM matrix */ 489 PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 490 PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 491 /* lambda_s^T F_UU w_1 */ 492 PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 493 PetscCall(VecResetArray(ts->vec_sensip_col)); 494 PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 495 /* lambda_s^T F_UP w_2 */ 496 PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 497 for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */ 498 PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 499 PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift)); 500 PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj])); 501 if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj])); 502 } 503 /* Solve stage equation LHS X = RHS for second-order adjoint */ 504 for (nadj = 0; nadj < ts->numcost; nadj++) { 505 KSPConvergedReason kspreason; 506 PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj])); 507 PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 508 if (kspreason < 0) { 509 ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 510 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj)); 511 } 512 } 513 } 514 515 /* Update sensitivities, and evaluate integrals if there is any */ 516 if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 517 th->shift = 1. / ((th->Theta - 1.) * adjoint_time_step); 518 th->stage_time = adjoint_ptime; 519 PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre)); 520 PetscCall(KSPSetOperators(ksp, J, Jpre)); 521 /* R_U at t_n */ 522 if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL)); 523 for (nadj = 0; nadj < ts->numcost; nadj++) { 524 PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj])); 525 if (quadJ) { 526 PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 527 PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 528 PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col)); 529 PetscCall(VecResetArray(ts->vec_drdu_col)); 530 PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 531 } 532 PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift)); 533 } 534 535 /* Second-order adjoint */ 536 if (ts->vecs_sensi2) { /* U_n */ 537 /* Get w1 at t_n from TLM matrix */ 538 PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr)); 539 PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 540 /* lambda_s^T F_UU w_1 */ 541 PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu)); 542 PetscCall(VecResetArray(ts->vec_sensip_col)); 543 PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr)); 544 /* lambda_s^T F_UU w_2 */ 545 PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup)); 546 for (nadj = 0; nadj < ts->numcost; nadj++) { 547 /* 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 */ 548 PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj])); 549 PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj])); 550 if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj])); 551 PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift)); 552 } 553 } 554 555 th->stage_time = ts->ptime; /* recover the old value */ 556 557 if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 558 /* U_{n+1} */ 559 th->shift = 1.0 / (adjoint_time_step * th->Theta); 560 PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol)); 561 PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE)); 562 if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 563 for (nadj = 0; nadj < ts->numcost; nadj++) { 564 PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 565 PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj])); 566 if (quadJp) { 567 PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 568 PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 569 PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col)); 570 PetscCall(VecResetArray(ts->vec_drdp_col)); 571 PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 572 } 573 } 574 if (ts->vecs_sensi2p) { /* second-order */ 575 /* Get w1 at t_{n+1} from TLM matrix */ 576 PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr)); 577 PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 578 /* lambda_s^T F_PU w_1 */ 579 PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 580 PetscCall(VecResetArray(ts->vec_sensip_col)); 581 PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 582 583 /* lambda_s^T F_PP w_2 */ 584 PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 585 for (nadj = 0; nadj < ts->numcost; nadj++) { 586 /* 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) */ 587 PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 588 PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj])); 589 if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj])); 590 if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj])); 591 } 592 } 593 594 /* U_s */ 595 PetscCall(VecZeroEntries(th->Xdot)); 596 PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE)); 597 if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp)); 598 for (nadj = 0; nadj < ts->numcost; nadj++) { 599 PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 600 PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj])); 601 if (quadJp) { 602 PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 603 PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 604 PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col)); 605 PetscCall(VecResetArray(ts->vec_drdp_col)); 606 PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 607 } 608 if (ts->vecs_sensi2p) { /* second-order */ 609 /* Get w1 at t_n from TLM matrix */ 610 PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr)); 611 PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 612 /* lambda_s^T F_PU w_1 */ 613 PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu)); 614 PetscCall(VecResetArray(ts->vec_sensip_col)); 615 PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr)); 616 /* lambda_s^T F_PP w_2 */ 617 PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp)); 618 for (nadj = 0; nadj < ts->numcost; nadj++) { 619 /* 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) */ 620 PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj])); 621 PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj])); 622 if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj])); 623 if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj])); 624 } 625 } 626 } 627 } 628 } else { /* one-stage case */ 629 th->shift = 0.0; 630 PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */ 631 PetscCall(KSPSetOperators(ksp, J, Jpre)); 632 if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 633 for (nadj = 0; nadj < ts->numcost; nadj++) { 634 PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj])); 635 PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj])); 636 if (quadJ) { 637 PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr)); 638 PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr)); 639 PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col)); 640 PetscCall(VecResetArray(ts->vec_drdu_col)); 641 PetscCall(MatDenseRestoreColumn(quadJ, &xarr)); 642 } 643 } 644 if (ts->vecs_sensip) { 645 th->shift = 1.0 / (adjoint_time_step * th->Theta); 646 PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); 647 PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 648 if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp)); 649 for (nadj = 0; nadj < ts->numcost; nadj++) { 650 PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj])); 651 PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj])); 652 if (quadJp) { 653 PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr)); 654 PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr)); 655 PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col)); 656 PetscCall(VecResetArray(ts->vec_drdp_col)); 657 PetscCall(MatDenseRestoreColumn(quadJp, &xarr)); 658 } 659 } 660 } 661 } 662 663 th->status = TS_STEP_COMPLETE; 664 PetscFunctionReturn(PETSC_SUCCESS); 665 } 666 667 static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X) 668 { 669 TS_Theta *th = (TS_Theta *)ts->data; 670 PetscReal dt = t - ts->ptime; 671 672 PetscFunctionBegin; 673 PetscCall(VecCopy(ts->vec_sol, th->X)); 674 if (th->endpoint) dt *= th->Theta; 675 PetscCall(VecWAXPY(X, dt, th->Xdot, th->X)); 676 PetscFunctionReturn(PETSC_SUCCESS); 677 } 678 679 static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) 680 { 681 TS_Theta *th = (TS_Theta *)ts->data; 682 Vec X = ts->vec_sol; /* X = solution */ 683 Vec Y = th->vec_lte_work; /* Y = X + LTE */ 684 PetscReal wltea, wlter; 685 686 PetscFunctionBegin; 687 if (!th->vec_sol_prev) { 688 *wlte = -1; 689 PetscFunctionReturn(PETSC_SUCCESS); 690 } 691 /* Cannot compute LTE in first step or in restart after event */ 692 if (ts->steprestart) { 693 *wlte = -1; 694 PetscFunctionReturn(PETSC_SUCCESS); 695 } 696 /* Compute LTE using backward differences with non-constant time step */ 697 { 698 PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 699 PetscReal a = 1 + h_prev / h; 700 PetscScalar scal[3]; 701 Vec vecs[3]; 702 scal[0] = +1 / a; 703 scal[1] = -1 / (a - 1); 704 scal[2] = +1 / (a * (a - 1)); 705 vecs[0] = X; 706 vecs[1] = th->X0; 707 vecs[2] = th->vec_sol_prev; 708 PetscCall(VecCopy(X, Y)); 709 PetscCall(VecMAXPY(Y, 3, scal, vecs)); 710 PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter)); 711 } 712 if (order) *order = 2; 713 PetscFunctionReturn(PETSC_SUCCESS); 714 } 715 716 static PetscErrorCode TSRollBack_Theta(TS ts) 717 { 718 TS_Theta *th = (TS_Theta *)ts->data; 719 TS quadts = ts->quadraturets; 720 721 PetscFunctionBegin; 722 if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol)); 723 th->status = TS_STEP_INCOMPLETE; 724 if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN)); 725 if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN)); 726 PetscFunctionReturn(PETSC_SUCCESS); 727 } 728 729 static PetscErrorCode TSForwardStep_Theta(TS ts) 730 { 731 TS_Theta *th = (TS_Theta *)ts->data; 732 TS quadts = ts->quadraturets; 733 Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 734 Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 735 PetscInt ntlm; 736 KSP ksp; 737 Mat J, Jpre, quadJ = NULL, quadJp = NULL; 738 PetscScalar *barr, *xarr; 739 PetscReal previous_shift; 740 741 PetscFunctionBegin; 742 previous_shift = th->shift; 743 PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN)); 744 745 if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN)); 746 PetscCall(SNESGetKSP(ts->snes, &ksp)); 747 PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL)); 748 if (quadts) { 749 PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL)); 750 PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL)); 751 } 752 753 /* Build RHS */ 754 if (th->endpoint) { /* 2-stage method*/ 755 th->shift = 1. / ((th->Theta - 1.) * th->time_step0); 756 PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 757 PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip)); 758 PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta)); 759 760 /* Add the f_p forcing terms */ 761 if (ts->Jacp) { 762 PetscCall(VecZeroEntries(th->Xdot)); 763 PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 764 PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN)); 765 th->shift = previous_shift; 766 PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol)); 767 PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 768 PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN)); 769 } 770 } else { /* 1-stage method */ 771 th->shift = 0.0; 772 PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 773 PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip)); 774 PetscCall(MatScale(MatDeltaFwdSensip, -1.)); 775 776 /* Add the f_p forcing terms */ 777 if (ts->Jacp) { 778 th->shift = previous_shift; 779 PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); 780 PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE)); 781 PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN)); 782 } 783 } 784 785 /* Build LHS */ 786 th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */ 787 if (th->endpoint) { 788 PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 789 } else { 790 PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE)); 791 } 792 PetscCall(KSPSetOperators(ksp, J, Jpre)); 793 794 /* 795 Evaluate the first stage of integral gradients with the 2-stage method: 796 drdu|t_n*S(t_n) + drdp|t_n 797 This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 798 */ 799 if (th->endpoint) { /* 2-stage method only */ 800 if (quadts && quadts->mat_sensip) { 801 PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL)); 802 PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp)); 803 PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp)); 804 PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 805 PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 806 } 807 } 808 809 /* Solve the tangent linear equation for forward sensitivities to parameters */ 810 for (ntlm = 0; ntlm < th->num_tlm; ntlm++) { 811 KSPConvergedReason kspreason; 812 PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr)); 813 PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr)); 814 if (th->endpoint) { 815 PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr)); 816 PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 817 PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col)); 818 PetscCall(VecResetArray(ts->vec_sensip_col)); 819 PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr)); 820 } else { 821 PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol)); 822 } 823 PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 824 if (kspreason < 0) { 825 ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 826 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm)); 827 } 828 PetscCall(VecResetArray(VecDeltaFwdSensipCol)); 829 PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr)); 830 } 831 832 /* 833 Evaluate the second stage of integral gradients with the 2-stage method: 834 drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 835 */ 836 if (quadts && quadts->mat_sensip) { 837 if (!th->endpoint) { 838 PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */ 839 PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL)); 840 PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp)); 841 PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp)); 842 PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 843 PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 844 PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); 845 } else { 846 PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL)); 847 PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp)); 848 PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp)); 849 PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN)); 850 PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN)); 851 } 852 } else { 853 if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); 854 } 855 PetscFunctionReturn(PETSC_SUCCESS); 856 } 857 858 static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[]) 859 { 860 TS_Theta *th = (TS_Theta *)ts->data; 861 862 PetscFunctionBegin; 863 if (ns) { 864 if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */ 865 else *ns = 2; /* endpoint form */ 866 } 867 if (stagesensip) { 868 if (!th->endpoint && th->Theta != 1.0) { 869 th->MatFwdStages[0] = th->MatDeltaFwdSensip; 870 } else { 871 th->MatFwdStages[0] = th->MatFwdSensip0; 872 th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */ 873 } 874 *stagesensip = th->MatFwdStages; 875 } 876 PetscFunctionReturn(PETSC_SUCCESS); 877 } 878 879 /*------------------------------------------------------------*/ 880 static PetscErrorCode TSReset_Theta(TS ts) 881 { 882 TS_Theta *th = (TS_Theta *)ts->data; 883 884 PetscFunctionBegin; 885 PetscCall(VecDestroy(&th->X)); 886 PetscCall(VecDestroy(&th->Xdot)); 887 PetscCall(VecDestroy(&th->X0)); 888 PetscCall(VecDestroy(&th->affine)); 889 890 PetscCall(VecDestroy(&th->vec_sol_prev)); 891 PetscCall(VecDestroy(&th->vec_lte_work)); 892 893 PetscCall(VecDestroy(&th->VecCostIntegral0)); 894 PetscFunctionReturn(PETSC_SUCCESS); 895 } 896 897 static PetscErrorCode TSAdjointReset_Theta(TS ts) 898 { 899 TS_Theta *th = (TS_Theta *)ts->data; 900 901 PetscFunctionBegin; 902 PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam)); 903 PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu)); 904 PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2)); 905 PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2)); 906 PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp)); 907 PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp)); 908 PetscFunctionReturn(PETSC_SUCCESS); 909 } 910 911 static PetscErrorCode TSDestroy_Theta(TS ts) 912 { 913 PetscFunctionBegin; 914 PetscCall(TSReset_Theta(ts)); 915 if (ts->dm) { 916 PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts)); 917 PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts)); 918 } 919 PetscCall(PetscFree(ts->data)); 920 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL)); 921 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL)); 922 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL)); 923 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL)); 924 PetscFunctionReturn(PETSC_SUCCESS); 925 } 926 927 /* 928 This defines the nonlinear equation that is to be solved with SNES 929 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 930 931 Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true, 932 otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is 933 U = (U_{n+1} + U0)/2 934 */ 935 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts) 936 { 937 TS_Theta *th = (TS_Theta *)ts->data; 938 Vec X0, Xdot; 939 DM dm, dmsave; 940 PetscReal shift = th->shift; 941 942 PetscFunctionBegin; 943 PetscCall(SNESGetDM(snes, &dm)); 944 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 945 PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot)); 946 if (x != X0) { 947 PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); 948 } else { 949 PetscCall(VecZeroEntries(Xdot)); 950 } 951 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 952 dmsave = ts->dm; 953 ts->dm = dm; 954 PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE)); 955 ts->dm = dmsave; 956 PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 957 PetscFunctionReturn(PETSC_SUCCESS); 958 } 959 960 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts) 961 { 962 TS_Theta *th = (TS_Theta *)ts->data; 963 Vec Xdot; 964 DM dm, dmsave; 965 PetscReal shift = th->shift; 966 967 PetscFunctionBegin; 968 PetscCall(SNESGetDM(snes, &dm)); 969 /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 970 PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot)); 971 972 dmsave = ts->dm; 973 ts->dm = dm; 974 PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE)); 975 ts->dm = dmsave; 976 PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 977 PetscFunctionReturn(PETSC_SUCCESS); 978 } 979 980 static PetscErrorCode TSForwardSetUp_Theta(TS ts) 981 { 982 TS_Theta *th = (TS_Theta *)ts->data; 983 TS quadts = ts->quadraturets; 984 985 PetscFunctionBegin; 986 /* combine sensitivities to parameters and sensitivities to initial values into one array */ 987 th->num_tlm = ts->num_parameters; 988 PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip)); 989 if (quadts && quadts->mat_sensip) { 990 PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp)); 991 PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0)); 992 } 993 /* backup sensitivity results for roll-backs */ 994 PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0)); 995 996 PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol)); 997 PetscFunctionReturn(PETSC_SUCCESS); 998 } 999 1000 static PetscErrorCode TSForwardReset_Theta(TS ts) 1001 { 1002 TS_Theta *th = (TS_Theta *)ts->data; 1003 TS quadts = ts->quadraturets; 1004 1005 PetscFunctionBegin; 1006 if (quadts && quadts->mat_sensip) { 1007 PetscCall(MatDestroy(&th->MatIntegralSensipTemp)); 1008 PetscCall(MatDestroy(&th->MatIntegralSensip0)); 1009 } 1010 PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol)); 1011 PetscCall(MatDestroy(&th->MatDeltaFwdSensip)); 1012 PetscCall(MatDestroy(&th->MatFwdSensip0)); 1013 PetscFunctionReturn(PETSC_SUCCESS); 1014 } 1015 1016 static PetscErrorCode TSSetUp_Theta(TS ts) 1017 { 1018 TS_Theta *th = (TS_Theta *)ts->data; 1019 TS quadts = ts->quadraturets; 1020 PetscBool match; 1021 1022 PetscFunctionBegin; 1023 if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */ 1024 PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0)); 1025 } 1026 if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X)); 1027 if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot)); 1028 if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0)); 1029 if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine)); 1030 1031 th->order = (th->Theta == 0.5) ? 2 : 1; 1032 th->shift = 1 / (th->Theta * ts->time_step); 1033 1034 PetscCall(TSGetDM(ts, &ts->dm)); 1035 PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts)); 1036 PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts)); 1037 1038 PetscCall(TSGetAdapt(ts, &ts->adapt)); 1039 PetscCall(TSAdaptCandidatesClear(ts->adapt)); 1040 PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match)); 1041 if (!match) { 1042 PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev)); 1043 PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work)); 1044 } 1045 PetscCall(TSGetSNES(ts, &ts->snes)); 1046 1047 ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE; 1048 PetscFunctionReturn(PETSC_SUCCESS); 1049 } 1050 1051 /*------------------------------------------------------------*/ 1052 1053 static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1054 { 1055 TS_Theta *th = (TS_Theta *)ts->data; 1056 1057 PetscFunctionBegin; 1058 PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam)); 1059 PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp)); 1060 if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu)); 1061 if (ts->vecs_sensi2) { 1062 PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2)); 1063 PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp)); 1064 /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 1065 if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu; 1066 if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup; 1067 } 1068 if (ts->vecs_sensi2p) { 1069 PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2)); 1070 /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 1071 if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu; 1072 if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp; 1073 } 1074 PetscFunctionReturn(PETSC_SUCCESS); 1075 } 1076 1077 static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems *PetscOptionsObject) 1078 { 1079 TS_Theta *th = (TS_Theta *)ts->data; 1080 1081 PetscFunctionBegin; 1082 PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options"); 1083 { 1084 PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL)); 1085 PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL)); 1086 PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL)); 1087 } 1088 PetscOptionsHeadEnd(); 1089 PetscFunctionReturn(PETSC_SUCCESS); 1090 } 1091 1092 static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer) 1093 { 1094 TS_Theta *th = (TS_Theta *)ts->data; 1095 PetscBool iascii; 1096 1097 PetscFunctionBegin; 1098 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1099 if (iascii) { 1100 PetscCall(PetscViewerASCIIPrintf(viewer, " Theta=%g\n", (double)th->Theta)); 1101 PetscCall(PetscViewerASCIIPrintf(viewer, " Extrapolation=%s\n", th->extrapolate ? "yes" : "no")); 1102 } 1103 PetscFunctionReturn(PETSC_SUCCESS); 1104 } 1105 1106 static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta) 1107 { 1108 TS_Theta *th = (TS_Theta *)ts->data; 1109 1110 PetscFunctionBegin; 1111 *theta = th->Theta; 1112 PetscFunctionReturn(PETSC_SUCCESS); 1113 } 1114 1115 static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta) 1116 { 1117 TS_Theta *th = (TS_Theta *)ts->data; 1118 1119 PetscFunctionBegin; 1120 PetscCheck(theta > 0 && theta <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Theta %g not in range (0,1]", (double)theta); 1121 th->Theta = theta; 1122 th->order = (th->Theta == 0.5) ? 2 : 1; 1123 PetscFunctionReturn(PETSC_SUCCESS); 1124 } 1125 1126 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint) 1127 { 1128 TS_Theta *th = (TS_Theta *)ts->data; 1129 1130 PetscFunctionBegin; 1131 *endpoint = th->endpoint; 1132 PetscFunctionReturn(PETSC_SUCCESS); 1133 } 1134 1135 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg) 1136 { 1137 TS_Theta *th = (TS_Theta *)ts->data; 1138 1139 PetscFunctionBegin; 1140 th->endpoint = flg; 1141 PetscFunctionReturn(PETSC_SUCCESS); 1142 } 1143 1144 #if defined(PETSC_HAVE_COMPLEX) 1145 static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) 1146 { 1147 PetscComplex z = xr + xi * PETSC_i, f; 1148 TS_Theta *th = (TS_Theta *)ts->data; 1149 const PetscReal one = 1.0; 1150 1151 PetscFunctionBegin; 1152 f = (one + (one - th->Theta) * z) / (one - th->Theta * z); 1153 *yr = PetscRealPartComplex(f); 1154 *yi = PetscImaginaryPartComplex(f); 1155 PetscFunctionReturn(PETSC_SUCCESS); 1156 } 1157 #endif 1158 1159 static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[]) 1160 { 1161 TS_Theta *th = (TS_Theta *)ts->data; 1162 1163 PetscFunctionBegin; 1164 if (ns) { 1165 if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */ 1166 else *ns = 2; /* endpoint form */ 1167 } 1168 if (Y) { 1169 if (!th->endpoint && th->Theta != 1.0) { 1170 th->Stages[0] = th->X; 1171 } else { 1172 th->Stages[0] = th->X0; 1173 th->Stages[1] = ts->vec_sol; /* stiffly accurate */ 1174 } 1175 *Y = th->Stages; 1176 } 1177 PetscFunctionReturn(PETSC_SUCCESS); 1178 } 1179 1180 /* ------------------------------------------------------------ */ 1181 /*MC 1182 TSTHETA - DAE solver using the implicit Theta method 1183 1184 Level: beginner 1185 1186 Options Database Keys: 1187 + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1188 . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 1189 - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 1190 1191 Notes: 1192 .vb 1193 -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 1194 -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 1195 -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 1196 .ve 1197 1198 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. 1199 1200 The midpoint variant is cast as a 1-stage implicit Runge-Kutta method. 1201 1202 .vb 1203 Theta | Theta 1204 ------------- 1205 | 1 1206 .ve 1207 1208 For the default Theta=0.5, this is also known as the implicit midpoint rule. 1209 1210 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1211 1212 .vb 1213 0 | 0 0 1214 1 | 1-Theta Theta 1215 ------------------- 1216 | 1-Theta Theta 1217 .ve 1218 1219 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1220 1221 To apply a diagonally implicit RK method to DAE, the stage formula 1222 1223 $ Y_i = X + h sum_j a_ij Y'_j 1224 1225 is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1226 1227 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()` 1228 M*/ 1229 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1230 { 1231 TS_Theta *th; 1232 1233 PetscFunctionBegin; 1234 ts->ops->reset = TSReset_Theta; 1235 ts->ops->adjointreset = TSAdjointReset_Theta; 1236 ts->ops->destroy = TSDestroy_Theta; 1237 ts->ops->view = TSView_Theta; 1238 ts->ops->setup = TSSetUp_Theta; 1239 ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1240 ts->ops->adjointreset = TSAdjointReset_Theta; 1241 ts->ops->step = TSStep_Theta; 1242 ts->ops->interpolate = TSInterpolate_Theta; 1243 ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 1244 ts->ops->rollback = TSRollBack_Theta; 1245 ts->ops->resizeregister = TSResizeRegister_Theta; 1246 ts->ops->setfromoptions = TSSetFromOptions_Theta; 1247 ts->ops->snesfunction = SNESTSFormFunction_Theta; 1248 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1249 #if defined(PETSC_HAVE_COMPLEX) 1250 ts->ops->linearstability = TSComputeLinearStability_Theta; 1251 #endif 1252 ts->ops->getstages = TSGetStages_Theta; 1253 ts->ops->adjointstep = TSAdjointStep_Theta; 1254 ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1255 ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 1256 ts->default_adapt_type = TSADAPTNONE; 1257 1258 ts->ops->forwardsetup = TSForwardSetUp_Theta; 1259 ts->ops->forwardreset = TSForwardReset_Theta; 1260 ts->ops->forwardstep = TSForwardStep_Theta; 1261 ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1262 1263 ts->usessnes = PETSC_TRUE; 1264 1265 PetscCall(PetscNew(&th)); 1266 ts->data = (void *)th; 1267 1268 th->VecsDeltaLam = NULL; 1269 th->VecsDeltaMu = NULL; 1270 th->VecsSensiTemp = NULL; 1271 th->VecsSensi2Temp = NULL; 1272 1273 th->extrapolate = PETSC_FALSE; 1274 th->Theta = 0.5; 1275 th->order = 2; 1276 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta)); 1277 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta)); 1278 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta)); 1279 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta)); 1280 PetscFunctionReturn(PETSC_SUCCESS); 1281 } 1282 1283 /*@ 1284 TSThetaGetTheta - Get the abscissa of the stage in (0,1] for `TSTHETA` 1285 1286 Not Collective 1287 1288 Input Parameter: 1289 . ts - timestepping context 1290 1291 Output Parameter: 1292 . theta - stage abscissa 1293 1294 Level: advanced 1295 1296 Note: 1297 Use of this function is normally only required to hack `TSTHETA` to use a modified integration scheme. 1298 1299 .seealso: [](ch_ts), `TSThetaSetTheta()`, `TSTHETA` 1300 @*/ 1301 PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta) 1302 { 1303 PetscFunctionBegin; 1304 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1305 PetscAssertPointer(theta, 2); 1306 PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta)); 1307 PetscFunctionReturn(PETSC_SUCCESS); 1308 } 1309 1310 /*@ 1311 TSThetaSetTheta - Set the abscissa of the stage in (0,1] for `TSTHETA` 1312 1313 Not Collective 1314 1315 Input Parameters: 1316 + ts - timestepping context 1317 - theta - stage abscissa 1318 1319 Options Database Key: 1320 . -ts_theta_theta <theta> - set theta 1321 1322 Level: intermediate 1323 1324 .seealso: [](ch_ts), `TSThetaGetTheta()`, `TSTHETA`, `TSCN` 1325 @*/ 1326 PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta) 1327 { 1328 PetscFunctionBegin; 1329 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1330 PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta)); 1331 PetscFunctionReturn(PETSC_SUCCESS); 1332 } 1333 1334 /*@ 1335 TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA` 1336 1337 Not Collective 1338 1339 Input Parameter: 1340 . ts - timestepping context 1341 1342 Output Parameter: 1343 . endpoint - `PETSC_TRUE` when using the endpoint variant 1344 1345 Level: advanced 1346 1347 .seealso: [](ch_ts), `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN` 1348 @*/ 1349 PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint) 1350 { 1351 PetscFunctionBegin; 1352 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1353 PetscAssertPointer(endpoint, 2); 1354 PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint)); 1355 PetscFunctionReturn(PETSC_SUCCESS); 1356 } 1357 1358 /*@ 1359 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA` 1360 1361 Not Collective 1362 1363 Input Parameters: 1364 + ts - timestepping context 1365 - flg - `PETSC_TRUE` to use the endpoint variant 1366 1367 Options Database Key: 1368 . -ts_theta_endpoint <flg> - use the endpoint variant 1369 1370 Level: intermediate 1371 1372 .seealso: [](ch_ts), `TSTHETA`, `TSCN` 1373 @*/ 1374 PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg) 1375 { 1376 PetscFunctionBegin; 1377 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1378 PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg)); 1379 PetscFunctionReturn(PETSC_SUCCESS); 1380 } 1381 1382 /* 1383 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1384 * The creation functions for these specializations are below. 1385 */ 1386 1387 static PetscErrorCode TSSetUp_BEuler(TS ts) 1388 { 1389 TS_Theta *th = (TS_Theta *)ts->data; 1390 1391 PetscFunctionBegin; 1392 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"); 1393 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"); 1394 PetscCall(TSSetUp_Theta(ts)); 1395 PetscFunctionReturn(PETSC_SUCCESS); 1396 } 1397 1398 static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer) 1399 { 1400 PetscFunctionBegin; 1401 PetscFunctionReturn(PETSC_SUCCESS); 1402 } 1403 1404 /*MC 1405 TSBEULER - ODE solver using the implicit backward Euler method 1406 1407 Level: beginner 1408 1409 Note: 1410 `TSBEULER` is equivalent to `TSTHETA` with Theta=1.0 or `-ts_type theta -ts_theta_theta 1.0` 1411 1412 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA` 1413 M*/ 1414 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1415 { 1416 PetscFunctionBegin; 1417 PetscCall(TSCreate_Theta(ts)); 1418 PetscCall(TSThetaSetTheta(ts, 1.0)); 1419 PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE)); 1420 ts->ops->setup = TSSetUp_BEuler; 1421 ts->ops->view = TSView_BEuler; 1422 PetscFunctionReturn(PETSC_SUCCESS); 1423 } 1424 1425 static PetscErrorCode TSSetUp_CN(TS ts) 1426 { 1427 TS_Theta *th = (TS_Theta *)ts->data; 1428 1429 PetscFunctionBegin; 1430 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"); 1431 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"); 1432 PetscCall(TSSetUp_Theta(ts)); 1433 PetscFunctionReturn(PETSC_SUCCESS); 1434 } 1435 1436 static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer) 1437 { 1438 PetscFunctionBegin; 1439 PetscFunctionReturn(PETSC_SUCCESS); 1440 } 1441 1442 /*MC 1443 TSCN - ODE solver using the implicit Crank-Nicolson method. 1444 1445 Level: beginner 1446 1447 Notes: 1448 `TSCN` is equivalent to `TSTHETA` with Theta=0.5 and the "endpoint" option set. I.e. 1449 .vb 1450 -ts_type theta 1451 -ts_theta_theta 0.5 1452 -ts_theta_endpoint 1453 .ve 1454 1455 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`, `TSType`, 1456 M*/ 1457 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1458 { 1459 PetscFunctionBegin; 1460 PetscCall(TSCreate_Theta(ts)); 1461 PetscCall(TSThetaSetTheta(ts, 0.5)); 1462 PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE)); 1463 ts->ops->setup = TSSetUp_CN; 1464 ts->ops->view = TSView_CN; 1465 PetscFunctionReturn(PETSC_SUCCESS); 1466 } 1467