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 ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 304 ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 305 ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 306 ierr = SNESSetFunction(ts->snes,pseudo->func,SNESTSFormFunction,ts);CHKERRQ(ierr); 307 /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will 308 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 309 context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets 310 the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and 311 snes->jacP should be the TS. */ 312 { 313 Mat A,B; 314 PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 315 void *ctx; 316 ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr); 317 ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr); 318 } 319 PetscFunctionReturn(0); 320 } 321 /*------------------------------------------------------------*/ 322 323 #undef __FUNCT__ 324 #define __FUNCT__ "TSPseudoMonitorDefault" 325 PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 326 { 327 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 328 PetscErrorCode ierr; 329 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 330 331 PetscFunctionBegin; 332 if (!ctx) { 333 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 334 } 335 ierr = PetscViewerASCIIMonitorPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 336 if (!ctx) { 337 ierr = PetscViewerASCIIMonitorDestroy(&viewer);CHKERRQ(ierr); 338 } 339 PetscFunctionReturn(0); 340 } 341 342 #undef __FUNCT__ 343 #define __FUNCT__ "TSSetFromOptions_Pseudo" 344 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 345 { 346 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 347 PetscErrorCode ierr; 348 PetscBool flg = PETSC_FALSE; 349 PetscViewerASCIIMonitor viewer; 350 351 PetscFunctionBegin; 352 ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 353 ierr = PetscOptionsBool("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 354 if (flg) { 355 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 356 ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 357 } 358 flg = PETSC_FALSE; 359 ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 360 if (flg) { 361 ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 362 } 363 ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 364 ierr = PetscOptionsTail();CHKERRQ(ierr); 365 PetscFunctionReturn(0); 366 } 367 368 #undef __FUNCT__ 369 #define __FUNCT__ "TSView_Pseudo" 370 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 371 { 372 PetscFunctionBegin; 373 PetscFunctionReturn(0); 374 } 375 376 /* ----------------------------------------------------------------------------- */ 377 #undef __FUNCT__ 378 #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 379 /*@C 380 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 381 last timestep. 382 383 Logically Collective on TS 384 385 Input Parameters: 386 + ts - timestep context 387 . dt - user-defined function to verify timestep 388 - ctx - [optional] user-defined context for private data 389 for the timestep verification routine (may be PETSC_NULL) 390 391 Level: advanced 392 393 Calling sequence of func: 394 . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 395 396 . update - latest solution vector 397 . ctx - [optional] timestep context 398 . newdt - the timestep to use for the next step 399 . flag - flag indicating whether the last time step was acceptable 400 401 Notes: 402 The routine set here will be called by TSPseudoVerifyTimeStep() 403 during the timestepping process. 404 405 .keywords: timestep, pseudo, set, verify 406 407 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 408 @*/ 409 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx) 410 { 411 PetscErrorCode ierr; 412 413 PetscFunctionBegin; 414 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 415 ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool *),void *),(ts,dt,ctx));CHKERRQ(ierr); 416 PetscFunctionReturn(0); 417 } 418 419 #undef __FUNCT__ 420 #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 421 /*@ 422 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 423 dt when using the TSPseudoDefaultTimeStep() routine. 424 425 Logically Collective on TS 426 427 Input Parameters: 428 + ts - the timestep context 429 - inc - the scaling factor >= 1.0 430 431 Options Database Key: 432 $ -ts_pseudo_increment <increment> 433 434 Level: advanced 435 436 .keywords: timestep, pseudo, set, increment 437 438 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 439 @*/ 440 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 441 { 442 PetscErrorCode ierr; 443 444 PetscFunctionBegin; 445 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 446 PetscValidLogicalCollectiveReal(ts,inc,2); 447 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 448 PetscFunctionReturn(0); 449 } 450 451 #undef __FUNCT__ 452 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 453 /*@ 454 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 455 is computed via the formula 456 $ dt = initial_dt*initial_fnorm/current_fnorm 457 rather than the default update, 458 $ dt = current_dt*previous_fnorm/current_fnorm. 459 460 Logically Collective on TS 461 462 Input Parameter: 463 . ts - the timestep context 464 465 Options Database Key: 466 $ -ts_pseudo_increment_dt_from_initial_dt 467 468 Level: advanced 469 470 .keywords: timestep, pseudo, set, increment 471 472 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 473 @*/ 474 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 475 { 476 PetscErrorCode ierr; 477 478 PetscFunctionBegin; 479 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 480 ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 481 PetscFunctionReturn(0); 482 } 483 484 485 #undef __FUNCT__ 486 #define __FUNCT__ "TSPseudoSetTimeStep" 487 /*@C 488 TSPseudoSetTimeStep - Sets the user-defined routine to be 489 called at each pseudo-timestep to update the timestep. 490 491 Logically Collective on TS 492 493 Input Parameters: 494 + ts - timestep context 495 . dt - function to compute timestep 496 - ctx - [optional] user-defined context for private data 497 required by the function (may be PETSC_NULL) 498 499 Level: intermediate 500 501 Calling sequence of func: 502 . func (TS ts,PetscReal *newdt,void *ctx); 503 504 . newdt - the newly computed timestep 505 . ctx - [optional] timestep context 506 507 Notes: 508 The routine set here will be called by TSPseudoComputeTimeStep() 509 during the timestepping process. 510 511 .keywords: timestep, pseudo, set 512 513 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 514 @*/ 515 PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 516 { 517 PetscErrorCode ierr; 518 519 PetscFunctionBegin; 520 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 521 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr); 522 PetscFunctionReturn(0); 523 } 524 525 /* ----------------------------------------------------------------------------- */ 526 527 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 528 EXTERN_C_BEGIN 529 #undef __FUNCT__ 530 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 531 PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 532 { 533 TS_Pseudo *pseudo; 534 535 PetscFunctionBegin; 536 pseudo = (TS_Pseudo*)ts->data; 537 pseudo->verify = dt; 538 pseudo->verifyctx = ctx; 539 PetscFunctionReturn(0); 540 } 541 EXTERN_C_END 542 543 EXTERN_C_BEGIN 544 #undef __FUNCT__ 545 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 546 PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 547 { 548 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 549 550 PetscFunctionBegin; 551 pseudo->dt_increment = inc; 552 PetscFunctionReturn(0); 553 } 554 EXTERN_C_END 555 556 EXTERN_C_BEGIN 557 #undef __FUNCT__ 558 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 559 PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 560 { 561 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 562 563 PetscFunctionBegin; 564 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 565 PetscFunctionReturn(0); 566 } 567 EXTERN_C_END 568 569 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 570 EXTERN_C_BEGIN 571 #undef __FUNCT__ 572 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 573 PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 574 { 575 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 576 577 PetscFunctionBegin; 578 pseudo->dt = dt; 579 pseudo->dtctx = ctx; 580 PetscFunctionReturn(0); 581 } 582 EXTERN_C_END 583 584 /* ----------------------------------------------------------------------------- */ 585 /*MC 586 TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 587 588 This method solves equations of the form 589 590 $ F(X,Xdot) = 0 591 592 for steady state using the iteration 593 594 $ [G'] S = -F(X,0) 595 $ X += S 596 597 where 598 599 $ G(Y) = F(Y,(Y-X)/dt) 600 601 This is linearly-implicit Euler with the residual always evaluated "at steady 602 state". See note below. 603 604 Options database keys: 605 + -ts_pseudo_increment <real> - ratio of increase dt 606 - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 607 608 Level: beginner 609 610 References: 611 Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 612 C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 613 614 Notes: 615 The residual computed by this method includes the transient term (Xdot is computed instead of 616 always being zero), but since the prediction from the last step is always the solution from the 617 last step, on the first Newton iteration we have 618 619 $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 620 621 Therefore, the linear system solved by the first Newton iteration is equivalent to the one 622 described above and in the papers. If the user chooses to perform multiple Newton iterations, the 623 algorithm is no longer the one described in the referenced papers. 624 625 .seealso: TSCreate(), TS, TSSetType() 626 627 M*/ 628 EXTERN_C_BEGIN 629 #undef __FUNCT__ 630 #define __FUNCT__ "TSCreate_Pseudo" 631 PetscErrorCode TSCreate_Pseudo(TS ts) 632 { 633 TS_Pseudo *pseudo; 634 PetscErrorCode ierr; 635 636 PetscFunctionBegin; 637 ts->ops->reset = TSReset_Pseudo; 638 ts->ops->destroy = TSDestroy_Pseudo; 639 ts->ops->view = TSView_Pseudo; 640 641 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 642 ts->ops->setup = TSSetUp_Pseudo; 643 ts->ops->step = TSStep_Pseudo; 644 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 645 ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 646 ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 647 648 /* create the required nonlinear solver context */ 649 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 650 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 651 652 ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 653 ts->data = (void*)pseudo; 654 655 pseudo->dt_increment = 1.1; 656 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 657 pseudo->dt = TSPseudoDefaultTimeStep; 658 659 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 660 "TSPseudoSetVerifyTimeStep_Pseudo", 661 TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 662 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 663 "TSPseudoSetTimeStepIncrement_Pseudo", 664 TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 665 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 666 "TSPseudoIncrementDtFromInitialDt_Pseudo", 667 TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 668 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 669 "TSPseudoSetTimeStep_Pseudo", 670 TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 671 PetscFunctionReturn(0); 672 } 673 EXTERN_C_END 674 675 #undef __FUNCT__ 676 #define __FUNCT__ "TSPseudoDefaultTimeStep" 677 /*@C 678 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 679 Use with TSPseudoSetTimeStep(). 680 681 Collective on TS 682 683 Input Parameters: 684 . ts - the timestep context 685 . dtctx - unused timestep context 686 687 Output Parameter: 688 . newdt - the timestep to use for the next step 689 690 Level: advanced 691 692 .keywords: timestep, pseudo, default 693 694 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 695 @*/ 696 PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 697 { 698 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 699 PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 700 PetscErrorCode ierr; 701 702 PetscFunctionBegin; 703 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 704 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 705 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 706 if (pseudo->initial_fnorm == 0.0) { 707 /* first time through so compute initial function norm */ 708 pseudo->initial_fnorm = pseudo->fnorm; 709 fnorm_previous = pseudo->fnorm; 710 } 711 if (pseudo->fnorm == 0.0) { 712 *newdt = 1.e12*inc*ts->time_step; 713 } else if (pseudo->increment_dt_from_initial_dt) { 714 *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 715 } else { 716 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 717 } 718 pseudo->fnorm_previous = pseudo->fnorm; 719 PetscFunctionReturn(0); 720 } 721