1 /* 2 Code for Timestepping with implicit backwards Euler. 3 */ 4 #include <petsc/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 PetscReal dt_max; /* maximum time step */ 24 PetscBool increment_dt_from_initial_dt; 25 } TS_Pseudo; 26 27 /* ------------------------------------------------------------------------------*/ 28 29 #undef __FUNCT__ 30 #define __FUNCT__ "TSPseudoComputeTimeStep" 31 /*@C 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: developer 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: TSPseudoTimeStepDefault(), TSPseudoSetTimeStep() 52 @*/ 53 PetscErrorCode 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__ "TSPseudoVerifyTimeStepDefault" 69 /*@C 70 TSPseudoVerifyTimeStepDefault - 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 TSPseudoVerifyTimeStepDefault(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(), TSPseudoVerifyTimeStepDefault() 125 @*/ 126 PetscErrorCode 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 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) pseudo->dt_initial = ts->time_step; 153 ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr); 154 next_time_step = ts->time_step; 155 ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr); 156 for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 157 ts->time_step = next_time_step; 158 ierr = TSPreStep(ts);CHKERRQ(ierr); 159 ierr = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr); 160 ierr = SNESSolve(ts->snes,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 ierr = TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));CHKERRQ(ierr); 165 ts->snes_its += its; ts->ksp_its += lits; 166 ierr = PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);CHKERRQ(ierr); 167 pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 168 ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr); 169 if (stepok) break; 170 } 171 if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 172 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 173 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); 174 PetscFunctionReturn(0); 175 } 176 if (reject >= ts->max_reject) { 177 ts->reason = TS_DIVERGED_STEP_REJECTED; 178 ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 179 PetscFunctionReturn(0); 180 } 181 ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr); 182 ts->ptime += 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 = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr); 213 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr); 214 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr); 215 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr); 216 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr); 217 PetscFunctionReturn(0); 218 } 219 220 /*------------------------------------------------------------*/ 221 222 #undef __FUNCT__ 223 #define __FUNCT__ "TSPseudoGetXdot" 224 /* 225 Compute Xdot = (X^{n+1}-X^n)/dt) = 0 226 */ 227 static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 228 { 229 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 230 const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 231 PetscScalar *xdot; 232 PetscErrorCode ierr; 233 PetscInt i,n; 234 235 PetscFunctionBegin; 236 ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 237 ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr); 238 ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 239 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 240 for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]); 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,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,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 = (PetscViewer) dummy; 323 324 PetscFunctionBegin; 325 if (!viewer) { 326 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);CHKERRQ(ierr); 327 } 328 if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 329 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 330 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 331 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 332 } 333 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 334 ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm);CHKERRQ(ierr); 335 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 336 PetscFunctionReturn(0); 337 } 338 339 #undef __FUNCT__ 340 #define __FUNCT__ "TSSetFromOptions_Pseudo" 341 static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptions *PetscOptionsObject,TS ts) 342 { 343 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 344 PetscErrorCode ierr; 345 PetscBool flg = PETSC_FALSE; 346 PetscViewer viewer; 347 348 PetscFunctionBegin; 349 ierr = PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");CHKERRQ(ierr); 350 ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,NULL);CHKERRQ(ierr); 351 if (flg) { 352 ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr); 353 ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 354 } 355 flg = PETSC_FALSE; 356 ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr); 357 if (flg) { 358 ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 359 } 360 ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);CHKERRQ(ierr); 361 ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);CHKERRQ(ierr); 362 ierr = PetscOptionsTail();CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 366 #undef __FUNCT__ 367 #define __FUNCT__ "TSView_Pseudo" 368 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 369 { 370 PetscErrorCode ierr; 371 372 PetscFunctionBegin; 373 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 374 PetscFunctionReturn(0); 375 } 376 377 /* ----------------------------------------------------------------------------- */ 378 #undef __FUNCT__ 379 #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 380 /*@C 381 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 382 last timestep. 383 384 Logically Collective on TS 385 386 Input Parameters: 387 + ts - timestep context 388 . dt - user-defined function to verify timestep 389 - ctx - [optional] user-defined context for private data 390 for the timestep verification routine (may be NULL) 391 392 Level: advanced 393 394 Calling sequence of func: 395 . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 396 397 . update - latest solution vector 398 . ctx - [optional] timestep context 399 . newdt - the timestep to use for the next step 400 . flag - flag indicating whether the last time step was acceptable 401 402 Notes: 403 The routine set here will be called by TSPseudoVerifyTimeStep() 404 during the timestepping process. 405 406 .keywords: timestep, pseudo, set, verify 407 408 .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep() 409 @*/ 410 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx) 411 { 412 PetscErrorCode ierr; 413 414 PetscFunctionBegin; 415 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 416 ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr); 417 PetscFunctionReturn(0); 418 } 419 420 #undef __FUNCT__ 421 #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 422 /*@ 423 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 424 dt when using the TSPseudoTimeStepDefault() routine. 425 426 Logically Collective on TS 427 428 Input Parameters: 429 + ts - the timestep context 430 - inc - the scaling factor >= 1.0 431 432 Options Database Key: 433 . -ts_pseudo_increment <increment> 434 435 Level: advanced 436 437 .keywords: timestep, pseudo, set, increment 438 439 .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 440 @*/ 441 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 442 { 443 PetscErrorCode ierr; 444 445 PetscFunctionBegin; 446 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 447 PetscValidLogicalCollectiveReal(ts,inc,2); 448 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 449 PetscFunctionReturn(0); 450 } 451 452 #undef __FUNCT__ 453 #define __FUNCT__ "TSPseudoSetMaxTimeStep" 454 /*@ 455 TSPseudoSetMaxTimeStep - Sets the maximum time step 456 when using the TSPseudoTimeStepDefault() routine. 457 458 Logically Collective on TS 459 460 Input Parameters: 461 + ts - the timestep context 462 - maxdt - the maximum time step, use a non-positive value to deactivate 463 464 Options Database Key: 465 . -ts_pseudo_max_dt <increment> 466 467 Level: advanced 468 469 .keywords: timestep, pseudo, set 470 471 .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 472 @*/ 473 PetscErrorCode TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt) 474 { 475 PetscErrorCode ierr; 476 477 PetscFunctionBegin; 478 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 479 PetscValidLogicalCollectiveReal(ts,maxdt,2); 480 ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 481 PetscFunctionReturn(0); 482 } 483 484 #undef __FUNCT__ 485 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 486 /*@ 487 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 488 is computed via the formula 489 $ dt = initial_dt*initial_fnorm/current_fnorm 490 rather than the default update, 491 $ dt = current_dt*previous_fnorm/current_fnorm. 492 493 Logically Collective on TS 494 495 Input Parameter: 496 . ts - the timestep context 497 498 Options Database Key: 499 . -ts_pseudo_increment_dt_from_initial_dt 500 501 Level: advanced 502 503 .keywords: timestep, pseudo, set, increment 504 505 .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 506 @*/ 507 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 508 { 509 PetscErrorCode ierr; 510 511 PetscFunctionBegin; 512 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 513 ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 514 PetscFunctionReturn(0); 515 } 516 517 518 #undef __FUNCT__ 519 #define __FUNCT__ "TSPseudoSetTimeStep" 520 /*@C 521 TSPseudoSetTimeStep - Sets the user-defined routine to be 522 called at each pseudo-timestep to update the timestep. 523 524 Logically Collective on TS 525 526 Input Parameters: 527 + ts - timestep context 528 . dt - function to compute timestep 529 - ctx - [optional] user-defined context for private data 530 required by the function (may be NULL) 531 532 Level: intermediate 533 534 Calling sequence of func: 535 . func (TS ts,PetscReal *newdt,void *ctx); 536 537 . newdt - the newly computed timestep 538 . ctx - [optional] timestep context 539 540 Notes: 541 The routine set here will be called by TSPseudoComputeTimeStep() 542 during the timestepping process. 543 If not set then TSPseudoTimeStepDefault() is automatically used 544 545 .keywords: timestep, pseudo, set 546 547 .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep() 548 @*/ 549 PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx) 550 { 551 PetscErrorCode ierr; 552 553 PetscFunctionBegin; 554 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 555 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr); 556 PetscFunctionReturn(0); 557 } 558 559 /* ----------------------------------------------------------------------------- */ 560 561 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 562 #undef __FUNCT__ 563 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 564 PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx) 565 { 566 TS_Pseudo *pseudo; 567 568 PetscFunctionBegin; 569 pseudo = (TS_Pseudo*)ts->data; 570 pseudo->verify = dt; 571 pseudo->verifyctx = ctx; 572 PetscFunctionReturn(0); 573 } 574 575 #undef __FUNCT__ 576 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 577 PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 578 { 579 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 580 581 PetscFunctionBegin; 582 pseudo->dt_increment = inc; 583 PetscFunctionReturn(0); 584 } 585 586 #undef __FUNCT__ 587 #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo" 588 PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt) 589 { 590 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 591 592 PetscFunctionBegin; 593 pseudo->dt_max = maxdt; 594 PetscFunctionReturn(0); 595 } 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 599 PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 600 { 601 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 602 603 PetscFunctionBegin; 604 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 605 PetscFunctionReturn(0); 606 } 607 608 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 609 #undef __FUNCT__ 610 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 611 PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx) 612 { 613 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 614 615 PetscFunctionBegin; 616 pseudo->dt = dt; 617 pseudo->dtctx = ctx; 618 PetscFunctionReturn(0); 619 } 620 621 /* ----------------------------------------------------------------------------- */ 622 /*MC 623 TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 624 625 This method solves equations of the form 626 627 $ F(X,Xdot) = 0 628 629 for steady state using the iteration 630 631 $ [G'] S = -F(X,0) 632 $ X += S 633 634 where 635 636 $ G(Y) = F(Y,(Y-X)/dt) 637 638 This is linearly-implicit Euler with the residual always evaluated "at steady 639 state". See note below. 640 641 Options database keys: 642 + -ts_pseudo_increment <real> - ratio of increase dt 643 - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 644 645 Level: beginner 646 647 References: 648 Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 649 C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 650 651 Notes: 652 The residual computed by this method includes the transient term (Xdot is computed instead of 653 always being zero), but since the prediction from the last step is always the solution from the 654 last step, on the first Newton iteration we have 655 656 $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 657 658 Therefore, the linear system solved by the first Newton iteration is equivalent to the one 659 described above and in the papers. If the user chooses to perform multiple Newton iterations, the 660 algorithm is no longer the one described in the referenced papers. 661 662 .seealso: TSCreate(), TS, TSSetType() 663 664 M*/ 665 #undef __FUNCT__ 666 #define __FUNCT__ "TSCreate_Pseudo" 667 PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 668 { 669 TS_Pseudo *pseudo; 670 PetscErrorCode ierr; 671 SNES snes; 672 SNESType stype; 673 674 PetscFunctionBegin; 675 ts->ops->reset = TSReset_Pseudo; 676 ts->ops->destroy = TSDestroy_Pseudo; 677 ts->ops->view = TSView_Pseudo; 678 679 ts->ops->setup = TSSetUp_Pseudo; 680 ts->ops->step = TSStep_Pseudo; 681 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 682 ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 683 ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 684 685 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 686 ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 687 if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 688 689 ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr); 690 ts->data = (void*)pseudo; 691 692 pseudo->dt_increment = 1.1; 693 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 694 pseudo->dt = TSPseudoTimeStepDefault; 695 pseudo->fnorm = -1; 696 697 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 698 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 699 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr); 700 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 701 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 702 PetscFunctionReturn(0); 703 } 704 705 #undef __FUNCT__ 706 #define __FUNCT__ "TSPseudoTimeStepDefault" 707 /*@C 708 TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. 709 Use with TSPseudoSetTimeStep(). 710 711 Collective on TS 712 713 Input Parameters: 714 . ts - the timestep context 715 . dtctx - unused timestep context 716 717 Output Parameter: 718 . newdt - the timestep to use for the next step 719 720 Level: advanced 721 722 .keywords: timestep, pseudo, default 723 724 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 725 @*/ 726 PetscErrorCode TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx) 727 { 728 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 729 PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 730 PetscErrorCode ierr; 731 732 PetscFunctionBegin; 733 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 734 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 735 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 736 if (pseudo->fnorm_initial == 0.0) { 737 /* first time through so compute initial function norm */ 738 pseudo->fnorm_initial = pseudo->fnorm; 739 fnorm_previous = pseudo->fnorm; 740 } 741 if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 742 else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 743 else *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 744 if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max); 745 pseudo->fnorm_previous = pseudo->fnorm; 746 PetscFunctionReturn(0); 747 } 748