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