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