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 = next_time_step; 183 ts->steps++; 184 PetscFunctionReturn(0); 185 } 186 187 /*------------------------------------------------------------*/ 188 #undef __FUNCT__ 189 #define __FUNCT__ "TSReset_Pseudo" 190 static PetscErrorCode TSReset_Pseudo(TS ts) 191 { 192 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 193 PetscErrorCode ierr; 194 195 PetscFunctionBegin; 196 ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr); 197 ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr); 198 ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr); 199 PetscFunctionReturn(0); 200 } 201 202 #undef __FUNCT__ 203 #define __FUNCT__ "TSDestroy_Pseudo" 204 static PetscErrorCode TSDestroy_Pseudo(TS ts) 205 { 206 PetscErrorCode ierr; 207 208 PetscFunctionBegin; 209 ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 210 ierr = PetscFree(ts->data);CHKERRQ(ierr); 211 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 212 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",PETSC_NULL);CHKERRQ(ierr); 213 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);CHKERRQ(ierr); 214 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 215 PetscFunctionReturn(0); 216 } 217 218 /*------------------------------------------------------------*/ 219 220 #undef __FUNCT__ 221 #define __FUNCT__ "TSPseudoGetXdot" 222 /* 223 Compute Xdot = (X^{n+1}-X^n)/dt) = 0 224 */ 225 static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 226 { 227 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 228 const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 229 PetscScalar *xdot; 230 PetscErrorCode ierr; 231 PetscInt i,n; 232 233 PetscFunctionBegin; 234 ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 235 ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr); 236 ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 237 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 238 for (i=0; i<n; i++) { 239 xdot[i] = mdt*(xnp1[i] - xn[i]); 240 } 241 ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 242 ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr); 243 ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 244 *Xdot = pseudo->xdot; 245 PetscFunctionReturn(0); 246 } 247 248 #undef __FUNCT__ 249 #define __FUNCT__ "SNESTSFormFunction_Pseudo" 250 /* 251 The transient residual is 252 253 F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 254 255 or for ODE, 256 257 (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 258 259 This is the function that must be evaluated for transient simulation and for 260 finite difference Jacobians. On the first Newton step, this algorithm uses 261 a guess of U^{n+1} = U^n in which case the transient term vanishes and the 262 residual is actually the steady state residual. Pseudotransient 263 continuation as described in the literature is a linearly implicit 264 algorithm, it only takes this one Newton step with the steady state 265 residual, and then advances to the next time step. 266 */ 267 static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 268 { 269 Vec Xdot; 270 PetscErrorCode ierr; 271 272 PetscFunctionBegin; 273 ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 274 ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 280 /* 281 This constructs the Jacobian needed for SNES. For DAE, this is 282 283 dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 284 285 and for ODE: 286 287 J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 288 */ 289 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 290 { 291 Vec Xdot; 292 PetscErrorCode ierr; 293 294 PetscFunctionBegin; 295 ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 296 ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr); 297 PetscFunctionReturn(0); 298 } 299 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "TSSetUp_Pseudo" 303 static PetscErrorCode TSSetUp_Pseudo(TS ts) 304 { 305 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 306 PetscErrorCode ierr; 307 308 PetscFunctionBegin; 309 ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 310 ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 311 ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 312 PetscFunctionReturn(0); 313 } 314 /*------------------------------------------------------------*/ 315 316 #undef __FUNCT__ 317 #define __FUNCT__ "TSPseudoMonitorDefault" 318 PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 319 { 320 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 321 PetscErrorCode ierr; 322 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 323 324 PetscFunctionBegin; 325 if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 326 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 327 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 328 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 329 } 330 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 331 ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 332 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 333 PetscFunctionReturn(0); 334 } 335 336 #undef __FUNCT__ 337 #define __FUNCT__ "TSSetFromOptions_Pseudo" 338 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 339 { 340 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 341 PetscErrorCode ierr; 342 PetscBool flg = PETSC_FALSE; 343 PetscViewer viewer; 344 345 PetscFunctionBegin; 346 ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 347 ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 348 if (flg) { 349 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr); 350 ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 351 } 352 flg = PETSC_FALSE; 353 ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 354 if (flg) { 355 ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 356 } 357 ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 358 359 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 360 ierr = PetscOptionsTail();CHKERRQ(ierr); 361 PetscFunctionReturn(0); 362 } 363 364 #undef __FUNCT__ 365 #define __FUNCT__ "TSView_Pseudo" 366 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 367 { 368 PetscErrorCode ierr; 369 370 PetscFunctionBegin; 371 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 372 PetscFunctionReturn(0); 373 } 374 375 /* ----------------------------------------------------------------------------- */ 376 #undef __FUNCT__ 377 #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 378 /*@C 379 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 380 last timestep. 381 382 Logically Collective on TS 383 384 Input Parameters: 385 + ts - timestep context 386 . dt - user-defined function to verify timestep 387 - ctx - [optional] user-defined context for private data 388 for the timestep verification routine (may be PETSC_NULL) 389 390 Level: advanced 391 392 Calling sequence of func: 393 . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 394 395 . update - latest solution vector 396 . ctx - [optional] timestep context 397 . newdt - the timestep to use for the next step 398 . flag - flag indicating whether the last time step was acceptable 399 400 Notes: 401 The routine set here will be called by TSPseudoVerifyTimeStep() 402 during the timestepping process. 403 404 .keywords: timestep, pseudo, set, verify 405 406 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 407 @*/ 408 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx) 409 { 410 PetscErrorCode ierr; 411 412 PetscFunctionBegin; 413 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 414 ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool *),void *),(ts,dt,ctx));CHKERRQ(ierr); 415 PetscFunctionReturn(0); 416 } 417 418 #undef __FUNCT__ 419 #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 420 /*@ 421 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 422 dt when using the TSPseudoDefaultTimeStep() routine. 423 424 Logically Collective on TS 425 426 Input Parameters: 427 + ts - the timestep context 428 - inc - the scaling factor >= 1.0 429 430 Options Database Key: 431 $ -ts_pseudo_increment <increment> 432 433 Level: advanced 434 435 .keywords: timestep, pseudo, set, increment 436 437 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 438 @*/ 439 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 440 { 441 PetscErrorCode ierr; 442 443 PetscFunctionBegin; 444 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 445 PetscValidLogicalCollectiveReal(ts,inc,2); 446 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 447 PetscFunctionReturn(0); 448 } 449 450 #undef __FUNCT__ 451 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 452 /*@ 453 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 454 is computed via the formula 455 $ dt = initial_dt*initial_fnorm/current_fnorm 456 rather than the default update, 457 $ dt = current_dt*previous_fnorm/current_fnorm. 458 459 Logically Collective on TS 460 461 Input Parameter: 462 . ts - the timestep context 463 464 Options Database Key: 465 $ -ts_pseudo_increment_dt_from_initial_dt 466 467 Level: advanced 468 469 .keywords: timestep, pseudo, set, increment 470 471 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 472 @*/ 473 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 474 { 475 PetscErrorCode ierr; 476 477 PetscFunctionBegin; 478 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 479 ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 480 PetscFunctionReturn(0); 481 } 482 483 484 #undef __FUNCT__ 485 #define __FUNCT__ "TSPseudoSetTimeStep" 486 /*@C 487 TSPseudoSetTimeStep - Sets the user-defined routine to be 488 called at each pseudo-timestep to update the timestep. 489 490 Logically Collective on TS 491 492 Input Parameters: 493 + ts - timestep context 494 . dt - function to compute timestep 495 - ctx - [optional] user-defined context for private data 496 required by the function (may be PETSC_NULL) 497 498 Level: intermediate 499 500 Calling sequence of func: 501 . func (TS ts,PetscReal *newdt,void *ctx); 502 503 . newdt - the newly computed timestep 504 . ctx - [optional] timestep context 505 506 Notes: 507 The routine set here will be called by TSPseudoComputeTimeStep() 508 during the timestepping process. 509 510 .keywords: timestep, pseudo, set 511 512 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 513 @*/ 514 PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 515 { 516 PetscErrorCode ierr; 517 518 PetscFunctionBegin; 519 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 520 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr); 521 PetscFunctionReturn(0); 522 } 523 524 /* ----------------------------------------------------------------------------- */ 525 526 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 527 EXTERN_C_BEGIN 528 #undef __FUNCT__ 529 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 530 PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 531 { 532 TS_Pseudo *pseudo; 533 534 PetscFunctionBegin; 535 pseudo = (TS_Pseudo*)ts->data; 536 pseudo->verify = dt; 537 pseudo->verifyctx = ctx; 538 PetscFunctionReturn(0); 539 } 540 EXTERN_C_END 541 542 EXTERN_C_BEGIN 543 #undef __FUNCT__ 544 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 545 PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 546 { 547 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 548 549 PetscFunctionBegin; 550 pseudo->dt_increment = inc; 551 PetscFunctionReturn(0); 552 } 553 EXTERN_C_END 554 555 EXTERN_C_BEGIN 556 #undef __FUNCT__ 557 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 558 PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 559 { 560 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 561 562 PetscFunctionBegin; 563 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 564 PetscFunctionReturn(0); 565 } 566 EXTERN_C_END 567 568 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 569 EXTERN_C_BEGIN 570 #undef __FUNCT__ 571 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 572 PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 573 { 574 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 575 576 PetscFunctionBegin; 577 pseudo->dt = dt; 578 pseudo->dtctx = ctx; 579 PetscFunctionReturn(0); 580 } 581 EXTERN_C_END 582 583 /* ----------------------------------------------------------------------------- */ 584 /*MC 585 TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 586 587 This method solves equations of the form 588 589 $ F(X,Xdot) = 0 590 591 for steady state using the iteration 592 593 $ [G'] S = -F(X,0) 594 $ X += S 595 596 where 597 598 $ G(Y) = F(Y,(Y-X)/dt) 599 600 This is linearly-implicit Euler with the residual always evaluated "at steady 601 state". See note below. 602 603 Options database keys: 604 + -ts_pseudo_increment <real> - ratio of increase dt 605 - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 606 607 Level: beginner 608 609 References: 610 Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 611 C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 612 613 Notes: 614 The residual computed by this method includes the transient term (Xdot is computed instead of 615 always being zero), but since the prediction from the last step is always the solution from the 616 last step, on the first Newton iteration we have 617 618 $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 619 620 Therefore, the linear system solved by the first Newton iteration is equivalent to the one 621 described above and in the papers. If the user chooses to perform multiple Newton iterations, the 622 algorithm is no longer the one described in the referenced papers. 623 624 .seealso: TSCreate(), TS, TSSetType() 625 626 M*/ 627 EXTERN_C_BEGIN 628 #undef __FUNCT__ 629 #define __FUNCT__ "TSCreate_Pseudo" 630 PetscErrorCode TSCreate_Pseudo(TS ts) 631 { 632 TS_Pseudo *pseudo; 633 PetscErrorCode ierr; 634 SNES snes; 635 const SNESType stype; 636 637 PetscFunctionBegin; 638 ts->ops->reset = TSReset_Pseudo; 639 ts->ops->destroy = TSDestroy_Pseudo; 640 ts->ops->view = TSView_Pseudo; 641 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 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 649 ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 650 if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);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 pseudo->fnorm = -1; 659 660 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 661 "TSPseudoSetVerifyTimeStep_Pseudo", 662 TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 663 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 664 "TSPseudoSetTimeStepIncrement_Pseudo", 665 TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 666 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 667 "TSPseudoIncrementDtFromInitialDt_Pseudo", 668 TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 669 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 670 "TSPseudoSetTimeStep_Pseudo", 671 TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 672 PetscFunctionReturn(0); 673 } 674 EXTERN_C_END 675 676 #undef __FUNCT__ 677 #define __FUNCT__ "TSPseudoDefaultTimeStep" 678 /*@C 679 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 680 Use with TSPseudoSetTimeStep(). 681 682 Collective on TS 683 684 Input Parameters: 685 . ts - the timestep context 686 . dtctx - unused timestep context 687 688 Output Parameter: 689 . newdt - the timestep to use for the next step 690 691 Level: advanced 692 693 .keywords: timestep, pseudo, default 694 695 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 696 @*/ 697 PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 698 { 699 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 700 PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 701 PetscErrorCode ierr; 702 703 PetscFunctionBegin; 704 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 705 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 706 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 707 if (pseudo->fnorm_initial == 0.0) { 708 /* first time through so compute initial function norm */ 709 pseudo->fnorm_initial = pseudo->fnorm; 710 fnorm_previous = pseudo->fnorm; 711 } 712 if (pseudo->fnorm == 0.0) { 713 *newdt = 1.e12*inc*ts->time_step; 714 } else if (pseudo->increment_dt_from_initial_dt) { 715 *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 716 } else { 717 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 718 } 719 pseudo->fnorm_previous = pseudo->fnorm; 720 PetscFunctionReturn(0); 721 } 722