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