1 /* 2 Code for Timestepping with implicit backwards Euler. 3 */ 4 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5 6 typedef struct { 7 Vec update; /* work vector where new solution is formed */ 8 Vec func; /* work vector where F(t[i],u[i]) is stored */ 9 Vec xdot; /* work vector for time derivative of state */ 10 11 /* information used for Pseudo-timestepping */ 12 13 PetscErrorCode (*dt)(TS, PetscReal *, void *); /* compute next timestep, and related context */ 14 void *dtctx; 15 PetscErrorCode (*verify)(TS, Vec, void *, PetscReal *, PetscBool *); /* verify previous timestep and related context */ 16 void *verifyctx; 17 18 PetscReal fnorm_initial, fnorm; /* original and current norm of F(u) */ 19 PetscReal fnorm_previous; 20 21 PetscReal dt_initial; /* initial time-step */ 22 PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 23 PetscReal dt_max; /* maximum time step */ 24 PetscBool increment_dt_from_initial_dt; 25 PetscReal fatol, frtol; 26 } TS_Pseudo; 27 28 /* ------------------------------------------------------------------------------*/ 29 30 /*@C 31 TSPseudoComputeTimeStep - Computes the next timestep for a currently running 32 pseudo-timestepping process. 33 34 Collective 35 36 Input Parameter: 37 . ts - timestep context 38 39 Output Parameter: 40 . dt - newly computed timestep 41 42 Level: developer 43 44 Note: 45 The routine to be called here to compute the timestep should be 46 set by calling `TSPseudoSetTimeStep()`. 47 48 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()` 49 @*/ 50 PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt) 51 { 52 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 53 54 PetscFunctionBegin; 55 PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 56 PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx)); 57 PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 58 PetscFunctionReturn(PETSC_SUCCESS); 59 } 60 61 /* ------------------------------------------------------------------------------*/ 62 /*@C 63 TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 64 65 Collective 66 67 Input Parameters: 68 + ts - the timestep context 69 . dtctx - unused timestep context 70 - update - latest solution vector 71 72 Output Parameters: 73 + newdt - the timestep to use for the next step 74 - flag - flag indicating whether the last time step was acceptable 75 76 Level: advanced 77 78 Note: 79 This routine always returns a flag of 1, indicating an acceptable 80 timestep. 81 82 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()` 83 @*/ 84 PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag) 85 { 86 PetscFunctionBegin; 87 *flag = PETSC_TRUE; 88 PetscFunctionReturn(PETSC_SUCCESS); 89 } 90 91 /*@ 92 TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 93 94 Collective 95 96 Input Parameters: 97 + ts - timestep context 98 - update - latest solution vector 99 100 Output Parameters: 101 + dt - newly computed timestep (if it had to shrink) 102 - flag - indicates if current timestep was ok 103 104 Level: advanced 105 106 Notes: 107 The routine to be called here to compute the timestep should be 108 set by calling `TSPseudoSetVerifyTimeStep()`. 109 110 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()` 111 @*/ 112 PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag) 113 { 114 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 115 116 PetscFunctionBegin; 117 *flag = PETSC_TRUE; 118 if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag)); 119 PetscFunctionReturn(PETSC_SUCCESS); 120 } 121 122 /* --------------------------------------------------------------------------------*/ 123 124 static PetscErrorCode TSStep_Pseudo(TS ts) 125 { 126 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 127 PetscInt nits, lits, reject; 128 PetscBool stepok; 129 PetscReal next_time_step = ts->time_step; 130 131 PetscFunctionBegin; 132 if (ts->steps == 0) pseudo->dt_initial = ts->time_step; 133 PetscCall(VecCopy(ts->vec_sol, pseudo->update)); 134 PetscCall(TSPseudoComputeTimeStep(ts, &next_time_step)); 135 for (reject = 0; reject < ts->max_reject; reject++, ts->reject++) { 136 ts->time_step = next_time_step; 137 PetscCall(TSPreStage(ts, ts->ptime + ts->time_step)); 138 PetscCall(SNESSolve(ts->snes, NULL, pseudo->update)); 139 PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 140 PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 141 ts->snes_its += nits; 142 ts->ksp_its += lits; 143 PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &(pseudo->update))); 144 PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, pseudo->update, &stepok)); 145 if (!stepok) { 146 next_time_step = ts->time_step; 147 continue; 148 } 149 pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 150 PetscCall(TSPseudoVerifyTimeStep(ts, pseudo->update, &next_time_step, &stepok)); 151 if (stepok) break; 152 } 153 if (reject >= ts->max_reject) { 154 ts->reason = TS_DIVERGED_STEP_REJECTED; 155 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, reject)); 156 PetscFunctionReturn(PETSC_SUCCESS); 157 } 158 159 PetscCall(VecCopy(pseudo->update, ts->vec_sol)); 160 ts->ptime += ts->time_step; 161 ts->time_step = next_time_step; 162 163 if (pseudo->fnorm < 0) { 164 PetscCall(VecZeroEntries(pseudo->xdot)); 165 PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 166 PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 167 } 168 if (pseudo->fnorm < pseudo->fatol) { 169 ts->reason = TS_CONVERGED_PSEUDO_FATOL; 170 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol)); 171 PetscFunctionReturn(PETSC_SUCCESS); 172 } 173 if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) { 174 ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 175 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g / fnorm_initial %g < frtol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->fnorm_initial, (double)pseudo->fatol)); 176 PetscFunctionReturn(PETSC_SUCCESS); 177 } 178 PetscFunctionReturn(PETSC_SUCCESS); 179 } 180 181 /*------------------------------------------------------------*/ 182 static PetscErrorCode TSReset_Pseudo(TS ts) 183 { 184 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 185 186 PetscFunctionBegin; 187 PetscCall(VecDestroy(&pseudo->update)); 188 PetscCall(VecDestroy(&pseudo->func)); 189 PetscCall(VecDestroy(&pseudo->xdot)); 190 PetscFunctionReturn(PETSC_SUCCESS); 191 } 192 193 static PetscErrorCode TSDestroy_Pseudo(TS ts) 194 { 195 PetscFunctionBegin; 196 PetscCall(TSReset_Pseudo(ts)); 197 PetscCall(PetscFree(ts->data)); 198 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL)); 199 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL)); 200 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL)); 201 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL)); 202 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL)); 203 PetscFunctionReturn(PETSC_SUCCESS); 204 } 205 206 /*------------------------------------------------------------*/ 207 208 /* 209 Compute Xdot = (X^{n+1}-X^n)/dt) = 0 210 */ 211 static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot) 212 { 213 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 214 const PetscScalar mdt = 1.0 / ts->time_step, *xnp1, *xn; 215 PetscScalar *xdot; 216 PetscInt i, n; 217 218 PetscFunctionBegin; 219 *Xdot = NULL; 220 PetscCall(VecGetArrayRead(ts->vec_sol, &xn)); 221 PetscCall(VecGetArrayRead(X, &xnp1)); 222 PetscCall(VecGetArray(pseudo->xdot, &xdot)); 223 PetscCall(VecGetLocalSize(X, &n)); 224 for (i = 0; i < n; i++) xdot[i] = mdt * (xnp1[i] - xn[i]); 225 PetscCall(VecRestoreArrayRead(ts->vec_sol, &xn)); 226 PetscCall(VecRestoreArrayRead(X, &xnp1)); 227 PetscCall(VecRestoreArray(pseudo->xdot, &xdot)); 228 *Xdot = pseudo->xdot; 229 PetscFunctionReturn(PETSC_SUCCESS); 230 } 231 232 /* 233 The transient residual is 234 235 F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 236 237 or for ODE, 238 239 (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 240 241 This is the function that must be evaluated for transient simulation and for 242 finite difference Jacobians. On the first Newton step, this algorithm uses 243 a guess of U^{n+1} = U^n in which case the transient term vanishes and the 244 residual is actually the steady state residual. Pseudotransient 245 continuation as described in the literature is a linearly implicit 246 algorithm, it only takes this one Newton step with the steady state 247 residual, and then advances to the next time step. 248 */ 249 static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts) 250 { 251 Vec Xdot; 252 253 PetscFunctionBegin; 254 PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 255 PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, Xdot, Y, PETSC_FALSE)); 256 PetscFunctionReturn(PETSC_SUCCESS); 257 } 258 259 /* 260 This constructs the Jacobian needed for SNES. For DAE, this is 261 262 dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 263 264 and for ODE: 265 266 J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 267 */ 268 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts) 269 { 270 Vec Xdot; 271 272 PetscFunctionBegin; 273 PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 274 PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE)); 275 PetscFunctionReturn(PETSC_SUCCESS); 276 } 277 278 static PetscErrorCode TSSetUp_Pseudo(TS ts) 279 { 280 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 281 282 PetscFunctionBegin; 283 PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update)); 284 PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func)); 285 PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot)); 286 PetscFunctionReturn(PETSC_SUCCESS); 287 } 288 /*------------------------------------------------------------*/ 289 290 static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy) 291 { 292 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 293 PetscViewer viewer = (PetscViewer)dummy; 294 295 PetscFunctionBegin; 296 if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 297 PetscCall(VecZeroEntries(pseudo->xdot)); 298 PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 299 PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 300 } 301 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 302 PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm)); 303 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 304 PetscFunctionReturn(PETSC_SUCCESS); 305 } 306 307 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems *PetscOptionsObject) 308 { 309 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 310 PetscBool flg = PETSC_FALSE; 311 PetscViewer viewer; 312 313 PetscFunctionBegin; 314 PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options"); 315 PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL)); 316 if (flg) { 317 PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer)); 318 PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 319 } 320 flg = pseudo->increment_dt_from_initial_dt; 321 PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL)); 322 pseudo->increment_dt_from_initial_dt = flg; 323 PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL)); 324 PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL)); 325 PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL)); 326 PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL)); 327 PetscOptionsHeadEnd(); 328 PetscFunctionReturn(PETSC_SUCCESS); 329 } 330 331 static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer) 332 { 333 PetscBool isascii; 334 335 PetscFunctionBegin; 336 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 337 if (isascii) { 338 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 339 PetscCall(PetscViewerASCIIPrintf(viewer, " frtol - relative tolerance in function value %g\n", (double)pseudo->frtol)); 340 PetscCall(PetscViewerASCIIPrintf(viewer, " fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol)); 341 PetscCall(PetscViewerASCIIPrintf(viewer, " dt_initial - initial timestep %g\n", (double)pseudo->dt_initial)); 342 PetscCall(PetscViewerASCIIPrintf(viewer, " dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment)); 343 PetscCall(PetscViewerASCIIPrintf(viewer, " dt_max - maximum time %g\n", (double)pseudo->dt_max)); 344 } 345 PetscFunctionReturn(PETSC_SUCCESS); 346 } 347 348 /* ----------------------------------------------------------------------------- */ 349 /*@C 350 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 351 last timestep. 352 353 Logically Collective 354 355 Input Parameters: 356 + ts - timestep context 357 . dt - user-defined function to verify timestep 358 - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`) 359 360 Calling sequence of `func`: 361 $ PetscErrorCode func(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag); 362 + ts - the time-step context 363 . update - latest solution vector 364 . ctx - [optional] timestep context 365 . newdt - the timestep to use for the next step 366 - flag - flag indicating whether the last time step was acceptable 367 368 Level: advanced 369 370 Note: 371 The routine set here will be called by `TSPseudoVerifyTimeStep()` 372 during the timestepping process. 373 374 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 375 @*/ 376 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS, Vec, void *, PetscReal *, PetscBool *), void *ctx) 377 { 378 PetscFunctionBegin; 379 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 380 PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode(*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx)); 381 PetscFunctionReturn(PETSC_SUCCESS); 382 } 383 384 /*@ 385 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 386 dt when using the TSPseudoTimeStepDefault() routine. 387 388 Logically Collective 389 390 Input Parameters: 391 + ts - the timestep context 392 - inc - the scaling factor >= 1.0 393 394 Options Database Key: 395 . -ts_pseudo_increment <increment> - set pseudo increment 396 397 Level: advanced 398 399 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 400 @*/ 401 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) 402 { 403 PetscFunctionBegin; 404 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 405 PetscValidLogicalCollectiveReal(ts, inc, 2); 406 PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc)); 407 PetscFunctionReturn(PETSC_SUCCESS); 408 } 409 410 /*@ 411 TSPseudoSetMaxTimeStep - Sets the maximum time step 412 when using the TSPseudoTimeStepDefault() routine. 413 414 Logically Collective 415 416 Input Parameters: 417 + ts - the timestep context 418 - maxdt - the maximum time step, use a non-positive value to deactivate 419 420 Options Database Key: 421 . -ts_pseudo_max_dt <increment> - set pseudo max dt 422 423 Level: advanced 424 425 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 426 @*/ 427 PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) 428 { 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 431 PetscValidLogicalCollectiveReal(ts, maxdt, 2); 432 PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 433 PetscFunctionReturn(PETSC_SUCCESS); 434 } 435 436 /*@ 437 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 438 is computed via the formula 439 $ dt = initial_dt*initial_fnorm/current_fnorm 440 rather than the default update, 441 $ dt = current_dt*previous_fnorm/current_fnorm. 442 443 Logically Collective 444 445 Input Parameter: 446 . ts - the timestep context 447 448 Options Database Key: 449 . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial dt to determine increment 450 451 Level: advanced 452 453 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 454 @*/ 455 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 456 { 457 PetscFunctionBegin; 458 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 459 PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts)); 460 PetscFunctionReturn(PETSC_SUCCESS); 461 } 462 463 /*@C 464 TSPseudoSetTimeStep - Sets the user-defined routine to be 465 called at each pseudo-timestep to update the timestep. 466 467 Logically Collective 468 469 Input Parameters: 470 + ts - timestep context 471 . dt - function to compute timestep 472 - ctx - [optional] user-defined context for private data required by the function (may be `NULL`) 473 474 Calling sequence of `dt`: 475 $ PetscErrorCode dt(TS ts, PetscReal *newdt, void *ctx); 476 + newdt - the newly computed timestep 477 - ctx - [optional] timestep context 478 479 Level: intermediate 480 481 Notes: 482 The routine set here will be called by `TSPseudoComputeTimeStep()` 483 during the timestepping process. 484 485 If not set then `TSPseudoTimeStepDefault()` is automatically used 486 487 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 488 @*/ 489 PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS, PetscReal *, void *), void *ctx) 490 { 491 PetscFunctionBegin; 492 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 493 PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode(*)(TS, PetscReal *, void *), void *), (ts, dt, ctx)); 494 PetscFunctionReturn(PETSC_SUCCESS); 495 } 496 497 /* ----------------------------------------------------------------------------- */ 498 499 typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/ 500 static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) 501 { 502 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 503 504 PetscFunctionBegin; 505 pseudo->verify = dt; 506 pseudo->verifyctx = ctx; 507 PetscFunctionReturn(PETSC_SUCCESS); 508 } 509 510 static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) 511 { 512 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 513 514 PetscFunctionBegin; 515 pseudo->dt_increment = inc; 516 PetscFunctionReturn(PETSC_SUCCESS); 517 } 518 519 static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) 520 { 521 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 522 523 PetscFunctionBegin; 524 pseudo->dt_max = maxdt; 525 PetscFunctionReturn(PETSC_SUCCESS); 526 } 527 528 static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 529 { 530 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 531 532 PetscFunctionBegin; 533 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 534 PetscFunctionReturn(PETSC_SUCCESS); 535 } 536 537 typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/ 538 static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) 539 { 540 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 541 542 PetscFunctionBegin; 543 pseudo->dt = dt; 544 pseudo->dtctx = ctx; 545 PetscFunctionReturn(PETSC_SUCCESS); 546 } 547 548 /* ----------------------------------------------------------------------------- */ 549 /*MC 550 TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 551 552 This method solves equations of the form 553 554 $ F(X,Xdot) = 0 555 556 for steady state using the iteration 557 558 .vb 559 [G'] S = -F(X,0) 560 X += S 561 .ve 562 563 where 564 565 $ G(Y) = F(Y,(Y-X)/dt) 566 567 This is linearly-implicit Euler with the residual always evaluated "at steady 568 state". See note below. 569 570 Options Database Keys: 571 + -ts_pseudo_increment <real> - ratio of increase dt 572 . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 573 . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol 574 - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol 575 576 Level: beginner 577 578 Notes: 579 The residual computed by this method includes the transient term (Xdot is computed instead of 580 always being zero), but since the prediction from the last step is always the solution from the 581 last step, on the first Newton iteration we have 582 583 $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 584 585 Therefore, the linear system solved by the first Newton iteration is equivalent to the one 586 described above and in the papers. If the user chooses to perform multiple Newton iterations, the 587 algorithm is no longer the one described in the referenced papers. 588 589 References: 590 + * - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003. 591 - * - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 592 593 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()` 594 M*/ 595 PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 596 { 597 TS_Pseudo *pseudo; 598 SNES snes; 599 SNESType stype; 600 601 PetscFunctionBegin; 602 ts->ops->reset = TSReset_Pseudo; 603 ts->ops->destroy = TSDestroy_Pseudo; 604 ts->ops->view = TSView_Pseudo; 605 ts->ops->setup = TSSetUp_Pseudo; 606 ts->ops->step = TSStep_Pseudo; 607 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 608 ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 609 ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 610 ts->default_adapt_type = TSADAPTNONE; 611 ts->usessnes = PETSC_TRUE; 612 613 PetscCall(TSGetSNES(ts, &snes)); 614 PetscCall(SNESGetType(snes, &stype)); 615 if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY)); 616 617 PetscCall(PetscNew(&pseudo)); 618 ts->data = (void *)pseudo; 619 620 pseudo->dt = TSPseudoTimeStepDefault; 621 pseudo->dtctx = NULL; 622 pseudo->dt_increment = 1.1; 623 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 624 pseudo->fnorm = -1; 625 pseudo->fnorm_initial = -1; 626 pseudo->fnorm_previous = -1; 627 #if defined(PETSC_USE_REAL_SINGLE) 628 pseudo->fatol = 1.e-25; 629 pseudo->frtol = 1.e-5; 630 #else 631 pseudo->fatol = 1.e-50; 632 pseudo->frtol = 1.e-12; 633 #endif 634 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo)); 635 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo)); 636 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo)); 637 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo)); 638 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo)); 639 PetscFunctionReturn(PETSC_SUCCESS); 640 } 641 642 /*@C 643 TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. Use with `TSPseudoSetTimeStep()`. 644 645 Collective 646 647 Input Parameters: 648 + ts - the timestep context 649 - dtctx - unused timestep context 650 651 Output Parameter: 652 . newdt - the timestep to use for the next step 653 654 Level: advanced 655 656 .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO` 657 @*/ 658 PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) 659 { 660 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 661 PetscReal inc = pseudo->dt_increment; 662 663 PetscFunctionBegin; 664 PetscCall(VecZeroEntries(pseudo->xdot)); 665 PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 666 PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 667 if (pseudo->fnorm_initial < 0) { 668 /* first time through so compute initial function norm */ 669 pseudo->fnorm_initial = pseudo->fnorm; 670 pseudo->fnorm_previous = pseudo->fnorm; 671 } 672 if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step; 673 else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm; 674 else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm; 675 if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max); 676 pseudo->fnorm_previous = pseudo->fnorm; 677 PetscFunctionReturn(PETSC_SUCCESS); 678 } 679