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