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