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