1 /* 2 Code for Timestepping with implicit backwards Euler. 3 */ 4 #include "src/ts/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 rhs; /* work vector for RHS; vec_sol/dt */ 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*,PetscTruth*); /* verify previous timestep and related context */ 16 void *verifyctx; 17 18 PetscReal initial_fnorm,fnorm; /* original and current norm of F(u) */ 19 PetscReal fnorm_previous; 20 21 PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 22 PetscTruth increment_dt_from_initial_dt; 23 } TS_Pseudo; 24 25 /* ------------------------------------------------------------------------------*/ 26 27 #undef __FUNCT__ 28 #define __FUNCT__ "TSPseudoComputeTimeStep" 29 /*@ 30 TSPseudoComputeTimeStep - Computes the next timestep for a currently running 31 pseudo-timestepping process. 32 33 Collective on TS 34 35 Input Parameter: 36 . ts - timestep context 37 38 Output Parameter: 39 . dt - newly computed timestep 40 41 Level: advanced 42 43 Notes: 44 The routine to be called here to compute the timestep should be 45 set by calling TSPseudoSetTimeStep(). 46 47 .keywords: timestep, pseudo, compute 48 49 .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 50 @*/ 51 PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 52 { 53 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54 PetscErrorCode ierr; 55 56 PetscFunctionBegin; 57 ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 58 ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 59 ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 60 PetscFunctionReturn(0); 61 } 62 63 64 /* ------------------------------------------------------------------------------*/ 65 #undef __FUNCT__ 66 #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 67 /*@C 68 TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 69 70 Collective on TS 71 72 Input Parameters: 73 + ts - the timestep context 74 . dtctx - unused timestep context 75 - update - latest solution vector 76 77 Output Parameters: 78 + newdt - the timestep to use for the next step 79 - flag - flag indicating whether the last time step was acceptable 80 81 Level: advanced 82 83 Note: 84 This routine always returns a flag of 1, indicating an acceptable 85 timestep. 86 87 .keywords: timestep, pseudo, default, verify 88 89 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 90 @*/ 91 PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscTruth *flag) 92 { 93 PetscFunctionBegin; 94 *flag = PETSC_TRUE; 95 PetscFunctionReturn(0); 96 } 97 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "TSPseudoVerifyTimeStep" 101 /*@ 102 TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 103 104 Collective on TS 105 106 Input Parameters: 107 + ts - timestep context 108 - update - latest solution vector 109 110 Output Parameters: 111 + dt - newly computed timestep (if it had to shrink) 112 - flag - indicates if current timestep was ok 113 114 Level: advanced 115 116 Notes: 117 The routine to be called here to compute the timestep should be 118 set by calling TSPseudoSetVerifyTimeStep(). 119 120 .keywords: timestep, pseudo, verify 121 122 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 123 @*/ 124 PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscTruth *flag) 125 { 126 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 127 PetscErrorCode ierr; 128 129 PetscFunctionBegin; 130 if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);} 131 132 ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 133 134 PetscFunctionReturn(0); 135 } 136 137 /* --------------------------------------------------------------------------------*/ 138 139 #undef __FUNCT__ 140 #define __FUNCT__ "TSStep_Pseudo" 141 static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime) 142 { 143 Vec sol = ts->vec_sol; 144 PetscErrorCode ierr; 145 PetscInt i,max_steps = ts->max_steps,its,lits; 146 PetscTruth ok; 147 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 148 PetscReal current_time_step; 149 150 PetscFunctionBegin; 151 *steps = -ts->steps; 152 153 ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr); 154 for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) { 155 ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 156 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 157 current_time_step = ts->time_step; 158 while (PETSC_TRUE) { 159 ts->ptime += current_time_step; 160 ierr = SNESSolve(ts->snes,pseudo->update);CHKERRQ(ierr); 161 ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr); 162 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 163 ts->nonlinear_its += its; ts->linear_its += lits; 164 ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr); 165 if (ok) break; 166 ts->ptime -= current_time_step; 167 current_time_step = ts->time_step; 168 } 169 ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr); 170 ts->steps++; 171 } 172 ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 173 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 174 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 175 176 *steps += ts->steps; 177 *ptime = ts->ptime; 178 PetscFunctionReturn(0); 179 } 180 181 /*------------------------------------------------------------*/ 182 #undef __FUNCT__ 183 #define __FUNCT__ "TSDestroy_Pseudo" 184 static PetscErrorCode TSDestroy_Pseudo(TS ts) 185 { 186 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 187 PetscErrorCode ierr; 188 189 PetscFunctionBegin; 190 if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);} 191 if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 192 if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 193 ierr = PetscFree(pseudo);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 198 /*------------------------------------------------------------*/ 199 200 /* 201 This defines the nonlinear equation that is to be solved with SNES 202 203 (U^{n+1} - U^{n})/dt - F(U^{n+1}) 204 */ 205 #undef __FUNCT__ 206 #define __FUNCT__ "TSPseudoFunction" 207 PetscErrorCode TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 208 { 209 TS ts = (TS) ctx; 210 PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 211 PetscErrorCode ierr; 212 PetscInt i,n; 213 214 PetscFunctionBegin; 215 /* apply user provided function */ 216 ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr); 217 /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 218 ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr); 219 ierr = VecGetArray(x,&unp1);CHKERRQ(ierr); 220 ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr); 221 ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 222 for (i=0; i<n; i++) { 223 Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 224 } 225 ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr); 226 ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr); 227 ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr); 228 229 PetscFunctionReturn(0); 230 } 231 232 /* 233 This constructs the Jacobian needed for SNES 234 235 J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 236 */ 237 #undef __FUNCT__ 238 #define __FUNCT__ "TSPseudoJacobian" 239 PetscErrorCode TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 240 { 241 TS ts = (TS) ctx; 242 PetscErrorCode ierr; 243 PetscScalar mone = -1.0,mdt = 1.0/ts->time_step; 244 245 PetscFunctionBegin; 246 /* construct users Jacobian */ 247 ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr); 248 249 /* shift and scale Jacobian */ 250 ierr = MatScale(&mone,*AA);CHKERRQ(ierr); 251 ierr = MatShift(&mdt,*AA);CHKERRQ(ierr); 252 if (*BB != *AA && *str != SAME_PRECONDITIONER) { 253 ierr = MatScale(&mone,*BB);CHKERRQ(ierr); 254 ierr = MatShift(&mdt,*BB);CHKERRQ(ierr); 255 } 256 257 PetscFunctionReturn(0); 258 } 259 260 261 #undef __FUNCT__ 262 #define __FUNCT__ "TSSetUp_Pseudo" 263 static PetscErrorCode TSSetUp_Pseudo(TS ts) 264 { 265 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 266 PetscErrorCode ierr; 267 268 PetscFunctionBegin; 269 /* ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); */ 270 ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 271 ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 272 ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 273 ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 274 PetscFunctionReturn(0); 275 } 276 /*------------------------------------------------------------*/ 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "TSPseudoDefaultMonitor" 280 PetscErrorCode TSPseudoDefaultMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 281 { 282 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 ierr = (*PetscHelpPrintf)(ts->comm,"TS %D dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 287 PetscFunctionReturn(0); 288 } 289 290 #undef __FUNCT__ 291 #define __FUNCT__ "TSSetFromOptions_Pseudo" 292 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 293 { 294 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 295 PetscErrorCode ierr; 296 PetscTruth flg; 297 298 PetscFunctionBegin; 299 300 ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 301 ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr); 302 if (flg) { 303 ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 304 } 305 ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr); 306 if (flg) { 307 ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 308 } 309 ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 310 ierr = PetscOptionsTail();CHKERRQ(ierr); 311 312 PetscFunctionReturn(0); 313 } 314 315 #undef __FUNCT__ 316 #define __FUNCT__ "TSView_Pseudo" 317 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 318 { 319 PetscFunctionBegin; 320 PetscFunctionReturn(0); 321 } 322 323 /* ----------------------------------------------------------------------------- */ 324 #undef __FUNCT__ 325 #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 326 /*@C 327 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 328 last timestep. 329 330 Collective on TS 331 332 Input Parameters: 333 + ts - timestep context 334 . dt - user-defined function to verify timestep 335 - ctx - [optional] user-defined context for private data 336 for the timestep verification routine (may be PETSC_NULL) 337 338 Level: advanced 339 340 Calling sequence of func: 341 . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscTruth *flag); 342 343 . update - latest solution vector 344 . ctx - [optional] timestep context 345 . newdt - the timestep to use for the next step 346 . flag - flag indicating whether the last time step was acceptable 347 348 Notes: 349 The routine set here will be called by TSPseudoVerifyTimeStep() 350 during the timestepping process. 351 352 .keywords: timestep, pseudo, set, verify 353 354 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 355 @*/ 356 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscTruth*),void* ctx) 357 { 358 PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscTruth *),void *); 359 360 PetscFunctionBegin; 361 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 362 363 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr); 364 if (f) { 365 ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 366 } 367 PetscFunctionReturn(0); 368 } 369 370 #undef __FUNCT__ 371 #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 372 /*@ 373 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 374 dt when using the TSPseudoDefaultTimeStep() routine. 375 376 Collective on TS 377 378 Input Parameters: 379 + ts - the timestep context 380 - inc - the scaling factor >= 1.0 381 382 Options Database Key: 383 $ -ts_pseudo_increment <increment> 384 385 Level: advanced 386 387 .keywords: timestep, pseudo, set, increment 388 389 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 390 @*/ 391 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 392 { 393 PetscErrorCode ierr,(*f)(TS,PetscReal); 394 395 PetscFunctionBegin; 396 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 397 398 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr); 399 if (f) { 400 ierr = (*f)(ts,inc);CHKERRQ(ierr); 401 } 402 PetscFunctionReturn(0); 403 } 404 405 #undef __FUNCT__ 406 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 407 /*@ 408 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 409 is computed via the formula 410 $ dt = initial_dt*initial_fnorm/current_fnorm 411 rather than the default update, 412 $ dt = current_dt*previous_fnorm/current_fnorm. 413 414 Collective on TS 415 416 Input Parameter: 417 . ts - the timestep context 418 419 Options Database Key: 420 $ -ts_pseudo_increment_dt_from_initial_dt 421 422 Level: advanced 423 424 .keywords: timestep, pseudo, set, increment 425 426 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 427 @*/ 428 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 429 { 430 PetscErrorCode ierr,(*f)(TS); 431 432 PetscFunctionBegin; 433 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 434 435 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr); 436 if (f) { 437 ierr = (*f)(ts);CHKERRQ(ierr); 438 } 439 PetscFunctionReturn(0); 440 } 441 442 443 #undef __FUNCT__ 444 #define __FUNCT__ "TSPseudoSetTimeStep" 445 /*@C 446 TSPseudoSetTimeStep - Sets the user-defined routine to be 447 called at each pseudo-timestep to update the timestep. 448 449 Collective on TS 450 451 Input Parameters: 452 + ts - timestep context 453 . dt - function to compute timestep 454 - ctx - [optional] user-defined context for private data 455 required by the function (may be PETSC_NULL) 456 457 Level: intermediate 458 459 Calling sequence of func: 460 . func (TS ts,PetscReal *newdt,void *ctx); 461 462 . newdt - the newly computed timestep 463 . ctx - [optional] timestep context 464 465 Notes: 466 The routine set here will be called by TSPseudoComputeTimeStep() 467 during the timestepping process. 468 469 .keywords: timestep, pseudo, set 470 471 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 472 @*/ 473 PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 474 { 475 PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *); 476 477 PetscFunctionBegin; 478 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 479 480 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr); 481 if (f) { 482 ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 483 } 484 PetscFunctionReturn(0); 485 } 486 487 /* ----------------------------------------------------------------------------- */ 488 489 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 490 EXTERN_C_BEGIN 491 #undef __FUNCT__ 492 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 493 PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 494 { 495 TS_Pseudo *pseudo; 496 497 PetscFunctionBegin; 498 pseudo = (TS_Pseudo*)ts->data; 499 pseudo->verify = dt; 500 pseudo->verifyctx = ctx; 501 PetscFunctionReturn(0); 502 } 503 EXTERN_C_END 504 505 EXTERN_C_BEGIN 506 #undef __FUNCT__ 507 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 508 PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 509 { 510 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 511 512 PetscFunctionBegin; 513 pseudo->dt_increment = inc; 514 PetscFunctionReturn(0); 515 } 516 EXTERN_C_END 517 518 EXTERN_C_BEGIN 519 #undef __FUNCT__ 520 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 521 PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 522 { 523 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 524 525 PetscFunctionBegin; 526 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 527 PetscFunctionReturn(0); 528 } 529 EXTERN_C_END 530 531 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 532 EXTERN_C_BEGIN 533 #undef __FUNCT__ 534 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 535 PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 536 { 537 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 538 539 PetscFunctionBegin; 540 pseudo->dt = dt; 541 pseudo->dtctx = ctx; 542 PetscFunctionReturn(0); 543 } 544 EXTERN_C_END 545 546 /* ----------------------------------------------------------------------------- */ 547 548 EXTERN_C_BEGIN 549 #undef __FUNCT__ 550 #define __FUNCT__ "TSCreate_Pseudo" 551 PetscErrorCode TSCreate_Pseudo(TS ts) 552 { 553 TS_Pseudo *pseudo; 554 PetscErrorCode ierr; 555 556 PetscFunctionBegin; 557 ts->ops->destroy = TSDestroy_Pseudo; 558 ts->ops->view = TSView_Pseudo; 559 560 if (ts->problem_type == TS_LINEAR) { 561 SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 562 } 563 if (!ts->A) { 564 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian"); 565 } 566 567 ts->ops->setup = TSSetUp_Pseudo; 568 ts->ops->step = TSStep_Pseudo; 569 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 570 571 /* create the required nonlinear solver context */ 572 ierr = SNESCreate(ts->comm,&ts->snes);CHKERRQ(ierr); 573 574 ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr); 575 PetscLogObjectMemory(ts,sizeof(TS_Pseudo)); 576 577 ierr = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr); 578 ts->data = (void*)pseudo; 579 580 pseudo->dt_increment = 1.1; 581 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 582 pseudo->dt = TSPseudoDefaultTimeStep; 583 584 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 585 "TSPseudoSetVerifyTimeStep_Pseudo", 586 TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 587 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 588 "TSPseudoSetTimeStepIncrement_Pseudo", 589 TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 590 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 591 "TSPseudoIncrementDtFromInitialDt_Pseudo", 592 TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 593 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 594 "TSPseudoSetTimeStep_Pseudo", 595 TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 596 PetscFunctionReturn(0); 597 } 598 EXTERN_C_END 599 600 #undef __FUNCT__ 601 #define __FUNCT__ "TSPseudoDefaultTimeStep" 602 /*@C 603 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 604 Use with TSPseudoSetTimeStep(). 605 606 Collective on TS 607 608 Input Parameters: 609 . ts - the timestep context 610 . dtctx - unused timestep context 611 612 Output Parameter: 613 . newdt - the timestep to use for the next step 614 615 Level: advanced 616 617 .keywords: timestep, pseudo, default 618 619 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 620 @*/ 621 PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 622 { 623 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 624 PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 625 PetscErrorCode ierr; 626 627 PetscFunctionBegin; 628 ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 629 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 630 if (pseudo->initial_fnorm == 0.0) { 631 /* first time through so compute initial function norm */ 632 pseudo->initial_fnorm = pseudo->fnorm; 633 fnorm_previous = pseudo->fnorm; 634 } 635 if (pseudo->fnorm == 0.0) { 636 *newdt = 1.e12*inc*ts->time_step; 637 } else if (pseudo->increment_dt_from_initial_dt) { 638 *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 639 } else { 640 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 641 } 642 pseudo->fnorm_previous = pseudo->fnorm; 643 PetscFunctionReturn(0); 644 } 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662