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