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