1 /* 2 Code for Timestepping with implicit backwards Euler. 3 */ 4 #include <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 initial_fnorm,fnorm; /* original and current norm of F(u) */ 19 PetscReal fnorm_previous; 20 21 PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 22 PetscBool increment_dt_from_initial_dt; 23 } TS_Pseudo; 24 25 /* ------------------------------------------------------------------------------*/ 26 27 #undef __FUNCT__ 28 #define __FUNCT__ "TSPseudoComputeTimeStep" 29 /*@ 30 TSPseudoComputeTimeStep - Computes the next timestep for a currently running 31 pseudo-timestepping process. 32 33 Collective on TS 34 35 Input Parameter: 36 . ts - timestep context 37 38 Output Parameter: 39 . dt - newly computed timestep 40 41 Level: advanced 42 43 Notes: 44 The routine to be called here to compute the timestep should be 45 set by calling TSPseudoSetTimeStep(). 46 47 .keywords: timestep, pseudo, compute 48 49 .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 50 @*/ 51 PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 52 { 53 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54 PetscErrorCode ierr; 55 56 PetscFunctionBegin; 57 ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 58 ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 59 ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 60 PetscFunctionReturn(0); 61 } 62 63 64 /* ------------------------------------------------------------------------------*/ 65 #undef __FUNCT__ 66 #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 67 /*@C 68 TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 69 70 Collective on TS 71 72 Input Parameters: 73 + ts - the timestep context 74 . dtctx - unused timestep context 75 - update - latest solution vector 76 77 Output Parameters: 78 + newdt - the timestep to use for the next step 79 - flag - flag indicating whether the last time step was acceptable 80 81 Level: advanced 82 83 Note: 84 This routine always returns a flag of 1, indicating an acceptable 85 timestep. 86 87 .keywords: timestep, pseudo, default, verify 88 89 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 90 @*/ 91 PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 92 { 93 PetscFunctionBegin; 94 *flag = PETSC_TRUE; 95 PetscFunctionReturn(0); 96 } 97 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "TSPseudoVerifyTimeStep" 101 /*@ 102 TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 103 104 Collective on TS 105 106 Input Parameters: 107 + ts - timestep context 108 - update - latest solution vector 109 110 Output Parameters: 111 + dt - newly computed timestep (if it had to shrink) 112 - flag - indicates if current timestep was ok 113 114 Level: advanced 115 116 Notes: 117 The routine to be called here to compute the timestep should be 118 set by calling TSPseudoSetVerifyTimeStep(). 119 120 .keywords: timestep, pseudo, verify 121 122 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 123 @*/ 124 PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 125 { 126 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 127 PetscErrorCode ierr; 128 129 PetscFunctionBegin; 130 if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);} 131 132 ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 133 134 PetscFunctionReturn(0); 135 } 136 137 /* --------------------------------------------------------------------------------*/ 138 139 #undef __FUNCT__ 140 #define __FUNCT__ "TSStep_Pseudo" 141 static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime) 142 { 143 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 144 Vec sol = ts->vec_sol, update = pseudo->update; 145 PetscErrorCode ierr; 146 PetscInt i,its,lits; 147 PetscBool ok; 148 PetscReal current_time_step; 149 150 PetscFunctionBegin; 151 *steps = -ts->steps; 152 *ptime = ts->ptime; 153 154 ierr = VecCopy(sol,update);CHKERRQ(ierr); 155 for (i=0; i<ts->max_steps && ts->ptime < ts->max_time; i++) { 156 157 ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 158 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 159 current_time_step = ts->time_step; 160 ierr = TSPreStep(ts);CHKERRQ(ierr); 161 while (PETSC_TRUE) { 162 ts->ptime += current_time_step; 163 ierr = SNESSolve(ts->snes,PETSC_NULL,update);CHKERRQ(ierr); 164 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 165 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 166 ts->nonlinear_its += its; ts->linear_its += lits; 167 ierr = TSPseudoVerifyTimeStep(ts,update,&ts->time_step,&ok);CHKERRQ(ierr); 168 if (ok) break; 169 ts->ptime -= current_time_step; 170 current_time_step = ts->time_step; 171 } 172 ierr = VecCopy(update,sol);CHKERRQ(ierr); 173 ts->steps++; 174 ierr = TSPostStep(ts);CHKERRQ(ierr); 175 } 176 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 177 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 178 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 179 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 180 181 *steps += ts->steps; 182 *ptime = ts->ptime; 183 PetscFunctionReturn(0); 184 } 185 186 /*------------------------------------------------------------*/ 187 #undef __FUNCT__ 188 #define __FUNCT__ "TSReset_Pseudo" 189 static PetscErrorCode TSReset_Pseudo(TS ts) 190 { 191 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 192 PetscErrorCode ierr; 193 194 PetscFunctionBegin; 195 ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr); 196 ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr); 197 ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr); 198 PetscFunctionReturn(0); 199 } 200 201 #undef __FUNCT__ 202 #define __FUNCT__ "TSDestroy_Pseudo" 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 PetscFunctionReturn(0); 211 } 212 213 /*------------------------------------------------------------*/ 214 215 #undef __FUNCT__ 216 #define __FUNCT__ "TSPseudoGetXdot" 217 /* 218 Compute Xdot = (X^{n+1}-X^n)/dt) = 0 219 */ 220 static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 221 { 222 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 223 PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn,*xdot; 224 PetscErrorCode ierr; 225 PetscInt i,n; 226 227 PetscFunctionBegin; 228 ierr = VecGetArray(ts->vec_sol,&xn);CHKERRQ(ierr); 229 ierr = VecGetArray(X,&xnp1);CHKERRQ(ierr); 230 ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 231 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 232 for (i=0; i<n; i++) { 233 xdot[i] = mdt*(xnp1[i] - xn[i]); 234 } 235 ierr = VecRestoreArray(ts->vec_sol,&xn);CHKERRQ(ierr); 236 ierr = VecRestoreArray(X,&xnp1);CHKERRQ(ierr); 237 ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 238 *Xdot = pseudo->xdot; 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "SNESTSFormFunction_Pseudo" 244 /* 245 The transient residual is 246 247 F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 248 249 or for ODE, 250 251 (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 252 253 This is the function that must be evaluated for transient simulation and for 254 finite difference Jacobians. On the first Newton step, this algorithm uses 255 a guess of U^{n+1} = U^n in which case the transient term vanishes and the 256 residual is actually the steady state residual. Pseudotransient 257 continuation as described in the literature is a linearly implicit 258 algorithm, it only takes this one Newton step with the steady state 259 residual, and then advances to the next time step. 260 */ 261 static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 262 { 263 Vec Xdot; 264 PetscErrorCode ierr; 265 266 PetscFunctionBegin; 267 ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 268 ierr = TSComputeIFunction(ts,ts->ptime,X,Xdot,Y);CHKERRQ(ierr); 269 PetscFunctionReturn(0); 270 } 271 272 #undef __FUNCT__ 273 #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 274 /* 275 This constructs the Jacobian needed for SNES. For DAE, this is 276 277 dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 278 279 and for ODE: 280 281 J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 282 */ 283 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 284 { 285 Vec Xdot; 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 290 ierr = TSComputeIJacobian(ts,ts->ptime,X,Xdot,1./ts->time_step,AA,BB,str);CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 295 #undef __FUNCT__ 296 #define __FUNCT__ "TSSetUp_Pseudo" 297 static PetscErrorCode TSSetUp_Pseudo(TS ts) 298 { 299 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 304 305 ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 306 ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 307 ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 308 ierr = SNESSetFunction(ts->snes,pseudo->func,SNESTSFormFunction,ts);CHKERRQ(ierr); 309 /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will 310 replace A and we don't want to mess with that. With -snes_mf, A and B will be replaced as well as the function and 311 context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets 312 the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and 313 snes->jacP should be the TS. */ 314 { 315 Mat A,B; 316 PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 317 void *ctx; 318 ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr); 319 ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr); 320 } 321 PetscFunctionReturn(0); 322 } 323 /*------------------------------------------------------------*/ 324 325 #undef __FUNCT__ 326 #define __FUNCT__ "TSPseudoMonitorDefault" 327 PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 328 { 329 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 330 PetscErrorCode ierr; 331 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 332 333 PetscFunctionBegin; 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","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 = PetscOptionsTail();CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 366 #undef __FUNCT__ 367 #define __FUNCT__ "TSView_Pseudo" 368 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 369 { 370 PetscFunctionBegin; 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 PETSC_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: TSPseudoDefaultVerifyTimeStep(), 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 TSPseudoDefaultTimeStep() 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(), TSPseudoDefaultTimeStep() 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__ "TSPseudoIncrementDtFromInitialDt" 451 /*@ 452 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 453 is computed via the formula 454 $ dt = initial_dt*initial_fnorm/current_fnorm 455 rather than the default update, 456 $ dt = current_dt*previous_fnorm/current_fnorm. 457 458 Logically Collective on TS 459 460 Input Parameter: 461 . ts - the timestep context 462 463 Options Database Key: 464 $ -ts_pseudo_increment_dt_from_initial_dt 465 466 Level: advanced 467 468 .keywords: timestep, pseudo, set, increment 469 470 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 471 @*/ 472 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 473 { 474 PetscErrorCode ierr; 475 476 PetscFunctionBegin; 477 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 478 ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 479 PetscFunctionReturn(0); 480 } 481 482 483 #undef __FUNCT__ 484 #define __FUNCT__ "TSPseudoSetTimeStep" 485 /*@C 486 TSPseudoSetTimeStep - Sets the user-defined routine to be 487 called at each pseudo-timestep to update the timestep. 488 489 Logically Collective on TS 490 491 Input Parameters: 492 + ts - timestep context 493 . dt - function to compute timestep 494 - ctx - [optional] user-defined context for private data 495 required by the function (may be PETSC_NULL) 496 497 Level: intermediate 498 499 Calling sequence of func: 500 . func (TS ts,PetscReal *newdt,void *ctx); 501 502 . newdt - the newly computed timestep 503 . ctx - [optional] timestep context 504 505 Notes: 506 The routine set here will be called by TSPseudoComputeTimeStep() 507 during the timestepping process. 508 509 .keywords: timestep, pseudo, set 510 511 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 512 @*/ 513 PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 514 { 515 PetscErrorCode ierr; 516 517 PetscFunctionBegin; 518 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 519 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr); 520 PetscFunctionReturn(0); 521 } 522 523 /* ----------------------------------------------------------------------------- */ 524 525 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 526 EXTERN_C_BEGIN 527 #undef __FUNCT__ 528 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 529 PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 530 { 531 TS_Pseudo *pseudo; 532 533 PetscFunctionBegin; 534 pseudo = (TS_Pseudo*)ts->data; 535 pseudo->verify = dt; 536 pseudo->verifyctx = ctx; 537 PetscFunctionReturn(0); 538 } 539 EXTERN_C_END 540 541 EXTERN_C_BEGIN 542 #undef __FUNCT__ 543 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 544 PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 545 { 546 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 547 548 PetscFunctionBegin; 549 pseudo->dt_increment = inc; 550 PetscFunctionReturn(0); 551 } 552 EXTERN_C_END 553 554 EXTERN_C_BEGIN 555 #undef __FUNCT__ 556 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 557 PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 558 { 559 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 560 561 PetscFunctionBegin; 562 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 563 PetscFunctionReturn(0); 564 } 565 EXTERN_C_END 566 567 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 568 EXTERN_C_BEGIN 569 #undef __FUNCT__ 570 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 571 PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 572 { 573 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 574 575 PetscFunctionBegin; 576 pseudo->dt = dt; 577 pseudo->dtctx = ctx; 578 PetscFunctionReturn(0); 579 } 580 EXTERN_C_END 581 582 /* ----------------------------------------------------------------------------- */ 583 /*MC 584 TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 585 586 This method solves equations of the form 587 588 $ F(X,Xdot) = 0 589 590 for steady state using the iteration 591 592 $ [G'] S = -F(X,0) 593 $ X += S 594 595 where 596 597 $ G(Y) = F(Y,(Y-X)/dt) 598 599 This is linearly-implicit Euler with the residual always evaluated "at steady 600 state". See note below. 601 602 Options database keys: 603 + -ts_pseudo_increment <real> - ratio of increase dt 604 - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 605 606 Level: beginner 607 608 References: 609 Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 610 C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 611 612 Notes: 613 The residual computed by this method includes the transient term (Xdot is computed instead of 614 always being zero), but since the prediction from the last step is always the solution from the 615 last step, on the first Newton iteration we have 616 617 $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 618 619 Therefore, the linear system solved by the first Newton iteration is equivalent to the one 620 described above and in the papers. If the user chooses to perform multiple Newton iterations, the 621 algorithm is no longer the one described in the referenced papers. 622 623 .seealso: TSCreate(), TS, TSSetType() 624 625 M*/ 626 EXTERN_C_BEGIN 627 #undef __FUNCT__ 628 #define __FUNCT__ "TSCreate_Pseudo" 629 PetscErrorCode TSCreate_Pseudo(TS ts) 630 { 631 TS_Pseudo *pseudo; 632 PetscErrorCode ierr; 633 634 PetscFunctionBegin; 635 ts->ops->reset = TSReset_Pseudo; 636 ts->ops->destroy = TSDestroy_Pseudo; 637 ts->ops->view = TSView_Pseudo; 638 639 ts->ops->setup = TSSetUp_Pseudo; 640 ts->ops->step = TSStep_Pseudo; 641 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 642 ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 643 ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 644 645 /* create the required nonlinear solver context */ 646 ts->problem_type = TS_NONLINEAR; 647 ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 648 649 ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 650 ts->data = (void*)pseudo; 651 652 pseudo->dt_increment = 1.1; 653 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 654 pseudo->dt = TSPseudoDefaultTimeStep; 655 656 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 657 "TSPseudoSetVerifyTimeStep_Pseudo", 658 TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 659 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 660 "TSPseudoSetTimeStepIncrement_Pseudo", 661 TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 662 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 663 "TSPseudoIncrementDtFromInitialDt_Pseudo", 664 TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 665 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 666 "TSPseudoSetTimeStep_Pseudo", 667 TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 668 PetscFunctionReturn(0); 669 } 670 EXTERN_C_END 671 672 #undef __FUNCT__ 673 #define __FUNCT__ "TSPseudoDefaultTimeStep" 674 /*@C 675 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 676 Use with TSPseudoSetTimeStep(). 677 678 Collective on TS 679 680 Input Parameters: 681 . ts - the timestep context 682 . dtctx - unused timestep context 683 684 Output Parameter: 685 . newdt - the timestep to use for the next step 686 687 Level: advanced 688 689 .keywords: timestep, pseudo, default 690 691 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 692 @*/ 693 PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 694 { 695 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 696 PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 697 PetscErrorCode ierr; 698 699 PetscFunctionBegin; 700 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 701 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 702 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 703 if (pseudo->initial_fnorm == 0.0) { 704 /* first time through so compute initial function norm */ 705 pseudo->initial_fnorm = pseudo->fnorm; 706 fnorm_previous = pseudo->fnorm; 707 } 708 if (pseudo->fnorm == 0.0) { 709 *newdt = 1.e12*inc*ts->time_step; 710 } else if (pseudo->increment_dt_from_initial_dt) { 711 *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 712 } else { 713 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 714 } 715 pseudo->fnorm_previous = pseudo->fnorm; 716 PetscFunctionReturn(0); 717 } 718