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