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