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