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