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