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 fnorm_initial,fnorm; /* original and current norm of F(u) */ 19 PetscReal fnorm_previous; 20 21 PetscReal dt_initial; /* initial time-step */ 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) 143 { 144 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 145 PetscInt its,lits,reject; 146 PetscBool stepok; 147 PetscReal next_time_step; 148 SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 149 PetscErrorCode ierr; 150 151 PetscFunctionBegin; 152 if (ts->steps == 0) { 153 pseudo->dt_initial = ts->time_step; 154 } 155 ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr); 156 next_time_step = ts->time_step; 157 ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr); 158 for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 159 ts->time_step = next_time_step; 160 ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->update);CHKERRQ(ierr); 161 ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 162 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 163 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 164 ts->nonlinear_its += its; ts->linear_its += lits; 165 ierr = PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);CHKERRQ(ierr); 166 pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 167 ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr); 168 if (stepok) break; 169 } 170 if (snesreason < 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 171 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 172 ierr = PetscInfo2(ts,"step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 173 PetscFunctionReturn(0); 174 } 175 if (reject >= ts->max_reject) { 176 ts->reason = TS_DIVERGED_STEP_REJECTED; 177 ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr); 181 ts->ptime += ts->time_step; 182 ts->time_step_prev = ts->time_step; 183 ts->time_step = next_time_step; 184 ts->steps++; 185 PetscFunctionReturn(0); 186 } 187 188 /*------------------------------------------------------------*/ 189 #undef __FUNCT__ 190 #define __FUNCT__ "TSReset_Pseudo" 191 static PetscErrorCode TSReset_Pseudo(TS ts) 192 { 193 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 194 PetscErrorCode ierr; 195 196 PetscFunctionBegin; 197 ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr); 198 ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr); 199 ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr); 200 PetscFunctionReturn(0); 201 } 202 203 #undef __FUNCT__ 204 #define __FUNCT__ "TSDestroy_Pseudo" 205 static PetscErrorCode TSDestroy_Pseudo(TS ts) 206 { 207 PetscErrorCode ierr; 208 209 PetscFunctionBegin; 210 ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 211 ierr = PetscFree(ts->data);CHKERRQ(ierr); 212 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 213 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",PETSC_NULL);CHKERRQ(ierr); 214 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);CHKERRQ(ierr); 215 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 216 PetscFunctionReturn(0); 217 } 218 219 /*------------------------------------------------------------*/ 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "TSPseudoGetXdot" 223 /* 224 Compute Xdot = (X^{n+1}-X^n)/dt) = 0 225 */ 226 static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 227 { 228 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 229 const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 230 PetscScalar *xdot; 231 PetscErrorCode ierr; 232 PetscInt i,n; 233 234 PetscFunctionBegin; 235 ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 236 ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr); 237 ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 238 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 239 for (i=0; i<n; i++) { 240 xdot[i] = mdt*(xnp1[i] - xn[i]); 241 } 242 ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 243 ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr); 244 ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 245 *Xdot = pseudo->xdot; 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "SNESTSFormFunction_Pseudo" 251 /* 252 The transient residual is 253 254 F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 255 256 or for ODE, 257 258 (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 259 260 This is the function that must be evaluated for transient simulation and for 261 finite difference Jacobians. On the first Newton step, this algorithm uses 262 a guess of U^{n+1} = U^n in which case the transient term vanishes and the 263 residual is actually the steady state residual. Pseudotransient 264 continuation as described in the literature is a linearly implicit 265 algorithm, it only takes this one Newton step with the steady state 266 residual, and then advances to the next time step. 267 */ 268 static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 269 { 270 Vec Xdot; 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 275 ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr); 276 PetscFunctionReturn(0); 277 } 278 279 #undef __FUNCT__ 280 #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 281 /* 282 This constructs the Jacobian needed for SNES. For DAE, this is 283 284 dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 285 286 and for ODE: 287 288 J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 289 */ 290 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 291 { 292 Vec Xdot; 293 PetscErrorCode ierr; 294 295 PetscFunctionBegin; 296 ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 297 ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr); 298 PetscFunctionReturn(0); 299 } 300 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "TSSetUp_Pseudo" 304 static PetscErrorCode TSSetUp_Pseudo(TS ts) 305 { 306 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 307 PetscErrorCode ierr; 308 309 PetscFunctionBegin; 310 ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 311 ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 312 ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 313 PetscFunctionReturn(0); 314 } 315 /*------------------------------------------------------------*/ 316 317 #undef __FUNCT__ 318 #define __FUNCT__ "TSPseudoMonitorDefault" 319 PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 320 { 321 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 322 PetscErrorCode ierr; 323 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 324 325 PetscFunctionBegin; 326 if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 327 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 328 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 329 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 330 } 331 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 332 ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 333 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 334 PetscFunctionReturn(0); 335 } 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "TSSetFromOptions_Pseudo" 339 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 340 { 341 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 342 PetscErrorCode ierr; 343 PetscBool flg = PETSC_FALSE; 344 PetscViewer viewer; 345 346 PetscFunctionBegin; 347 ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 348 ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 349 if (flg) { 350 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr); 351 ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 352 } 353 flg = PETSC_FALSE; 354 ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 355 if (flg) { 356 ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 357 } 358 ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 359 360 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 361 ierr = PetscOptionsTail();CHKERRQ(ierr); 362 PetscFunctionReturn(0); 363 } 364 365 #undef __FUNCT__ 366 #define __FUNCT__ "TSView_Pseudo" 367 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 368 { 369 PetscErrorCode ierr; 370 371 PetscFunctionBegin; 372 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 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 SNES snes; 636 const SNESType stype; 637 638 PetscFunctionBegin; 639 ts->ops->reset = TSReset_Pseudo; 640 ts->ops->destroy = TSDestroy_Pseudo; 641 ts->ops->view = TSView_Pseudo; 642 643 ts->ops->setup = TSSetUp_Pseudo; 644 ts->ops->step = TSStep_Pseudo; 645 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 646 ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 647 ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 648 649 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 650 ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 651 if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 652 653 ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 654 ts->data = (void*)pseudo; 655 656 pseudo->dt_increment = 1.1; 657 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 658 pseudo->dt = TSPseudoDefaultTimeStep; 659 pseudo->fnorm = -1; 660 661 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 662 "TSPseudoSetVerifyTimeStep_Pseudo", 663 TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 664 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 665 "TSPseudoSetTimeStepIncrement_Pseudo", 666 TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 667 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 668 "TSPseudoIncrementDtFromInitialDt_Pseudo", 669 TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 670 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 671 "TSPseudoSetTimeStep_Pseudo", 672 TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 673 PetscFunctionReturn(0); 674 } 675 EXTERN_C_END 676 677 #undef __FUNCT__ 678 #define __FUNCT__ "TSPseudoDefaultTimeStep" 679 /*@C 680 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 681 Use with TSPseudoSetTimeStep(). 682 683 Collective on TS 684 685 Input Parameters: 686 . ts - the timestep context 687 . dtctx - unused timestep context 688 689 Output Parameter: 690 . newdt - the timestep to use for the next step 691 692 Level: advanced 693 694 .keywords: timestep, pseudo, default 695 696 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 697 @*/ 698 PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 699 { 700 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 701 PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 702 PetscErrorCode ierr; 703 704 PetscFunctionBegin; 705 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 706 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 707 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 708 if (pseudo->fnorm_initial == 0.0) { 709 /* first time through so compute initial function norm */ 710 pseudo->fnorm_initial = pseudo->fnorm; 711 fnorm_previous = pseudo->fnorm; 712 } 713 if (pseudo->fnorm == 0.0) { 714 *newdt = 1.e12*inc*ts->time_step; 715 } else if (pseudo->increment_dt_from_initial_dt) { 716 *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 717 } else { 718 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 719 } 720 pseudo->fnorm_previous = pseudo->fnorm; 721 PetscFunctionReturn(0); 722 } 723