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