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