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