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