1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: posindep.c,v 1.25 1998/03/20 22:51:32 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(TS ts ) 174 { 175 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 176 int ierr; 177 178 PetscFunctionBegin; 179 ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 180 if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 181 if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 182 if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 183 PetscFree(pseudo); 184 PetscFunctionReturn(0); 185 } 186 187 188 /*------------------------------------------------------------*/ 189 /* 190 This matrix shell multiply where user provided Shell matrix 191 */ 192 193 #undef __FUNC__ 194 #define __FUNC__ "TSPseudoMatMult" 195 int TSPseudoMatMult(Mat mat,Vec x,Vec y) 196 { 197 TS ts; 198 Scalar mdt,mone = -1.0; 199 int ierr; 200 201 PetscFunctionBegin; 202 ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr); 203 mdt = 1.0/ts->time_step; 204 205 /* apply user provided function */ 206 ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 207 /* shift and scale by 1/dt - F */ 208 ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 /* 213 This defines the nonlinear equation that is to be solved with SNES 214 215 (U^{n+1} - U^{n})/dt - F(U^{n+1}) 216 */ 217 #undef __FUNC__ 218 #define __FUNC__ "TSPseudoFunction" 219 int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 220 { 221 TS ts = (TS) ctx; 222 Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 223 int ierr,i,n; 224 225 PetscFunctionBegin; 226 /* apply user provided function */ 227 ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 228 /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 229 ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 230 ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 231 ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 232 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 233 for ( i=0; i<n; i++ ) { 234 Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 235 } 236 ierr = VecRestoreArray(ts->vec_sol,&un); 237 ierr = VecRestoreArray(x,&unp1); 238 ierr = VecRestoreArray(y,&Funp1); 239 240 PetscFunctionReturn(0); 241 } 242 243 /* 244 This constructs the Jacobian needed for SNES 245 246 J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 247 */ 248 #undef __FUNC__ 249 #define __FUNC__ "TSPseudoJacobian" 250 int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 251 { 252 TS ts = (TS) ctx; 253 int ierr; 254 Scalar mone = -1.0, mdt = 1.0/ts->time_step; 255 MatType mtype; 256 257 PetscFunctionBegin; 258 /* construct users Jacobian */ 259 if (ts->rhsjacobian) { 260 ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 261 } 262 263 /* shift and scale Jacobian, if not a shell matrix */ 264 ierr = MatGetType(*AA,&mtype,PETSC_NULL); 265 if (mtype != MATSHELL) { 266 ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 267 ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 268 } 269 ierr = MatGetType(*BB,&mtype,PETSC_NULL); 270 if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 271 ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 272 ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 273 } 274 275 PetscFunctionReturn(0); 276 } 277 278 279 #undef __FUNC__ 280 #define __FUNC__ "TSSetUp_Pseudo" 281 static int TSSetUp_Pseudo(TS ts) 282 { 283 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 284 int ierr, M, m; 285 286 PetscFunctionBegin; 287 ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 288 ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 289 ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 290 if (ts->Ashell) { /* construct new shell matrix */ 291 ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 292 ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 293 ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 294 ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 295 } 296 ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 297 PetscFunctionReturn(0); 298 } 299 /*------------------------------------------------------------*/ 300 301 #undef __FUNC__ 302 #define __FUNC__ "TSPseudoDefaultMonitor" 303 int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 304 { 305 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 306 307 PetscFunctionBegin; 308 (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 309 PetscFunctionReturn(0); 310 } 311 312 #undef __FUNC__ 313 #define __FUNC__ "TSSetFromOptions_Pseudo" 314 static int TSSetFromOptions_Pseudo(TS ts) 315 { 316 int ierr,flg; 317 double inc; 318 319 PetscFunctionBegin; 320 ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 321 322 ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 323 if (flg) { 324 ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 325 } 326 ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 327 if (flg) { 328 ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 329 } 330 ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr); 331 if (flg) { 332 ierr = TSPseudoIncrementDtFromInitialDt(ts); CHKERRQ(ierr); 333 } 334 PetscFunctionReturn(0); 335 } 336 337 #undef __FUNC__ 338 #define __FUNC__ "TSPrintHelp_Pseudo" 339 static int TSPrintHelp_Pseudo(TS ts,char *p) 340 { 341 PetscFunctionBegin; 342 (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n"); 343 (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p); 344 (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p); 345 (*PetscHelpPrintf)(ts->comm," initial fnorm/current fnorm to determine new timestep\n"); 346 (*PetscHelpPrintf)(ts->comm," default is current_dt * previous fnorm/current fnorm\n"); 347 PetscFunctionReturn(0); 348 } 349 350 #undef __FUNC__ 351 #define __FUNC__ "TSView_Pseudo" 352 static int TSView_Pseudo(TS ts,Viewer viewer) 353 { 354 PetscFunctionBegin; 355 PetscFunctionReturn(0); 356 } 357 358 /* ----------------------------------------------------------------------------- */ 359 #undef __FUNC__ 360 #define __FUNC__ "TSPseudoSetVerifyTimeStep" 361 /*@ 362 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 363 last timestep. 364 365 Input Parameters: 366 . ts - timestep context 367 . dt - user-defined function to verify timestep 368 . ctx - [optional] user-defined context for private data 369 for the timestep verification routine (may be PETSC_NULL) 370 371 Calling sequence of func: 372 . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 373 374 . update - latest solution vector 375 . ctx - [optional] timestep context 376 . newdt - the timestep to use for the next step 377 . flag - flag indicating whether the last time step was acceptable 378 379 Notes: 380 The routine set here will be called by TSPseudoVerifyTimeStep() 381 during the timestepping process. 382 383 .keywords: timestep, pseudo, set, verify 384 385 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 386 @*/ 387 int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 388 { 389 int ierr, (*f)(TS,int (*)(TS,Vec,void*,double *,int *),void *); 390 391 PetscFunctionBegin; 392 PetscValidHeaderSpecific(ts,TS_COOKIE); 393 394 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep",(void **)&f);CHKERRQ(ierr); 395 if (f) { 396 ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 397 } 398 PetscFunctionReturn(0); 399 } 400 401 #undef __FUNC__ 402 #define __FUNC__ "TSPseudoSetTimeStepIncrement" 403 /*@ 404 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 405 dt when using the TSPseudoDefaultTimeStep() routine. 406 407 Input Parameters: 408 . ts - the timestep context 409 . inc - the scaling factor >= 1.0 410 411 Options Database Key: 412 $ -ts_pseudo_increment <increment> 413 414 .keywords: timestep, pseudo, set, increment 415 416 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 417 @*/ 418 int TSPseudoSetTimeStepIncrement(TS ts,double inc) 419 { 420 int ierr, (*f)(TS,double); 421 422 PetscFunctionBegin; 423 PetscValidHeaderSpecific(ts,TS_COOKIE); 424 425 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncremen",(void **)&f);CHKERRQ(ierr); 426 if (f) { 427 ierr = (*f)(ts,inc);CHKERRQ(ierr); 428 } 429 PetscFunctionReturn(0); 430 } 431 432 #undef __FUNC__ 433 #define __FUNC__ "TSPseudoIncrementDtFromInitialDt" 434 /*@ 435 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 436 is computed via the formula 437 $ dt = initial_dt*initial_fnorm/current_fnorm 438 rather than the default update, 439 $ dt = current_dt*previous_fnorm/current_fnorm. 440 441 Input Parameter: 442 . ts - the timestep context 443 444 Options Database Key: 445 $ -ts_pseudo_increment_dt_from_initial_dt 446 447 .keywords: timestep, pseudo, set, increment 448 449 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 450 @*/ 451 int TSPseudoIncrementDtFromInitialDt(TS ts) 452 { 453 int ierr, (*f)(TS); 454 455 PetscFunctionBegin; 456 PetscValidHeaderSpecific(ts,TS_COOKIE); 457 458 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt",(void **)&f);CHKERRQ(ierr); 459 if (f) { 460 ierr = (*f)(ts);CHKERRQ(ierr); 461 } 462 PetscFunctionReturn(0); 463 } 464 465 466 #undef __FUNC__ 467 #define __FUNC__ "TSPseudoSetTimeStep" 468 /*@ 469 TSPseudoSetTimeStep - Sets the user-defined routine to be 470 called at each pseudo-timestep to update the timestep. 471 472 Input Parameters: 473 . ts - timestep context 474 . dt - function to compute timestep 475 . ctx - [optional] user-defined context for private data 476 required by the function (may be PETSC_NULL) 477 478 Calling sequence of func: 479 . func (TS ts,double *newdt,void *ctx); 480 481 . newdt - the newly computed timestep 482 . ctx - [optional] timestep context 483 484 Notes: 485 The routine set here will be called by TSPseudoComputeTimeStep() 486 during the timestepping process. 487 488 .keywords: timestep, pseudo, set 489 490 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 491 @*/ 492 int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 493 { 494 int ierr, (*f)(TS,int (*)(TS,double *,void *),void *); 495 496 PetscFunctionBegin; 497 PetscValidHeaderSpecific(ts,TS_COOKIE); 498 499 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep",(void **)&f);CHKERRQ(ierr); 500 if (f) { 501 ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 502 } 503 PetscFunctionReturn(0); 504 } 505 506 /* ----------------------------------------------------------------------------- */ 507 508 #undef __FUNC__ 509 #define __FUNC__ "TSPseudoSetVerifyTimeStep_Pseudo" 510 int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 511 { 512 TS_Pseudo *pseudo; 513 514 PetscFunctionBegin; 515 pseudo = (TS_Pseudo*) ts->data; 516 pseudo->verify = dt; 517 pseudo->verifyctx = ctx; 518 PetscFunctionReturn(0); 519 } 520 521 #undef __FUNC__ 522 #define __FUNC__ "TSPseudoSetTimeStepIncrement_Pseudo" 523 int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,double inc) 524 { 525 TS_Pseudo *pseudo; 526 527 PetscFunctionBegin; 528 pseudo = (TS_Pseudo*) ts->data; 529 pseudo->dt_increment = inc; 530 PetscFunctionReturn(0); 531 } 532 533 #undef __FUNC__ 534 #define __FUNC__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 535 int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 536 { 537 TS_Pseudo *pseudo; 538 539 PetscFunctionBegin; 540 pseudo = (TS_Pseudo*) ts->data; 541 pseudo->increment_dt_from_initial_dt = 1; 542 PetscFunctionReturn(0); 543 } 544 545 #undef __FUNC__ 546 #define __FUNC__ "TSPseudoSetTimeStep_Pseudo" 547 int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,double*,void*),void* ctx) 548 { 549 TS_Pseudo *pseudo; 550 551 PetscFunctionBegin; 552 pseudo = (TS_Pseudo*) ts->data; 553 pseudo->dt = dt; 554 pseudo->dtctx = ctx; 555 PetscFunctionReturn(0); 556 } 557 558 /* ----------------------------------------------------------------------------- */ 559 560 #undef __FUNC__ 561 #define __FUNC__ "TSCreate_Pseudo" 562 int TSCreate_Pseudo(TS ts ) 563 { 564 TS_Pseudo *pseudo; 565 int ierr; 566 MatType mtype; 567 568 PetscFunctionBegin; 569 ts->destroy = TSDestroy_Pseudo; 570 ts->printhelp = TSPrintHelp_Pseudo; 571 ts->view = TSView_Pseudo; 572 573 if (ts->problem_type == TS_LINEAR) { 574 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems"); 575 } 576 if (!ts->A) { 577 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian"); 578 } 579 ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 580 if (mtype == MATSHELL) { 581 ts->Ashell = ts->A; 582 } 583 ts->setup = TSSetUp_Pseudo; 584 ts->step = TSStep_Pseudo; 585 ts->setfromoptions = TSSetFromOptions_Pseudo; 586 587 /* create the required nonlinear solver context */ 588 ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 589 590 pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 591 PLogObjectMemory(ts,sizeof(TS_Pseudo)); 592 593 PetscMemzero(pseudo,sizeof(TS_Pseudo)); 594 ts->data = (void *) pseudo; 595 596 pseudo->dt_increment = 1.1; 597 pseudo->increment_dt_from_initial_dt = 0; 598 pseudo->dt = TSPseudoDefaultTimeStep; 599 600 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep", 601 "TSPseudoSetVerifyTimeStep_Pseudo", 602 (void*)TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 603 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement", 604 "TSPseudoSetTimeStepIncrement_Pseudo", 605 (void*)TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 606 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt", 607 "TSPseudoIncrementDtFromInitialDt_Pseudo", 608 (void*)TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 609 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep","TSPseudoSetTimeStep_Pseudo", 610 (void*)TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 611 PetscFunctionReturn(0); 612 } 613 614 #undef __FUNC__ 615 #define __FUNC__ "TSPseudoDefaultTimeStep" 616 /*@C 617 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 618 Use with TSPseudoSetTimeStep(). 619 620 Input Parameters: 621 . ts - the timestep context 622 . dtctx - unused timestep context 623 624 Output Parameter: 625 . newdt - the timestep to use for the next step 626 627 .keywords: timestep, pseudo, default 628 629 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 630 @*/ 631 int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 632 { 633 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 634 double inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 635 int ierr; 636 637 PetscFunctionBegin; 638 ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 639 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 640 if (pseudo->initial_fnorm == 0.0) { 641 /* first time through so compute initial function norm */ 642 pseudo->initial_fnorm = pseudo->fnorm; 643 fnorm_previous = pseudo->fnorm; 644 } 645 if (pseudo->fnorm == 0.0) { 646 *newdt = 1.e12*inc*ts->time_step; 647 } 648 else if (pseudo->increment_dt_from_initial_dt) { 649 *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 650 } else { 651 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 652 } 653 pseudo->fnorm_previous = pseudo->fnorm; 654 PetscFunctionReturn(0); 655 } 656 657