1 #ifdef PETSC_RCS_HEADER 2 3 #endif 4 5 #include "src/ts/tsimpl.h" /*I "ts.h" I*/ 6 7 #undef __FUNC__ 8 #define __FUNC__ "TSComputeRHSFunction" 9 /* 10 TSComputeRHSFunction - Evaluates the right-hand-side function. 11 12 Note: If the user did not provide a function but merely a matrix, 13 this routine applies the matrix. 14 */ 15 int TSComputeRHSFunction(TS ts,double t,Vec x, Vec y) 16 { 17 int ierr; 18 19 PetscFunctionBegin; 20 PetscValidHeaderSpecific(ts,TS_COOKIE); 21 PetscValidHeader(x); PetscValidHeader(y); 22 23 if (ts->rhsfunction) { 24 PetscStackPush("TS user right-hand-side function"); 25 ierr = (*ts->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 26 PetscStackPop; 27 PetscFunctionReturn(0); 28 } 29 30 if (ts->rhsmatrix) { /* assemble matrix for this timestep */ 31 MatStructure flg; 32 PetscStackPush("TS user right-hand-side matrix function"); 33 ierr = (*ts->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 34 PetscStackPop; 35 } 36 ierr = MatMult(ts->A,x,y);CHKERRQ(ierr); 37 38 /* apply user-provided boundary conditions (only needed if these are time dependent) */ 39 ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr); 40 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNC__ 45 #define __FUNC__ "TSSetRHSFunction" 46 /*@C 47 TSSetRHSFunction - Sets the routine for evaluating the function, 48 F(t,u), where U_t = F(t,u). 49 50 Collective on TS 51 52 Input Parameters: 53 + ts - the TS context obtained from TSCreate() 54 . f - routine for evaluating the right-hand-side function 55 - ctx - [optional] user-defined context for private data for the 56 function evaluation routine (may be PETSC_NULL) 57 58 Calling sequence of func: 59 $ func (TS ts,double t,Vec u,Vec F,void *ctx); 60 61 + t - current timestep 62 . u - input vector 63 . F - function vector 64 - ctx - [optional] user-defined function context 65 66 Important: 67 The user MUST call either this routine or TSSetRHSMatrix(). 68 69 Level: beginner 70 71 .keywords: TS, timestep, set, right-hand-side, function 72 73 .seealso: TSSetRHSMatrix() 74 @*/ 75 int TSSetRHSFunction(TS ts,int (*f)(TS,double,Vec,Vec,void*),void *ctx) 76 { 77 PetscFunctionBegin; 78 79 PetscValidHeaderSpecific(ts,TS_COOKIE); 80 if (ts->problem_type == TS_LINEAR) { 81 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Cannot set function for linear problem"); 82 } 83 ts->rhsfunction = f; 84 ts->funP = ctx; 85 PetscFunctionReturn(0); 86 } 87 88 #undef __FUNC__ 89 #define __FUNC__ "TSSetRHSMatrix" 90 /*@C 91 TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U. 92 Also sets the location to store A. 93 94 Collective on TS 95 96 Input Parameters: 97 + ts - the TS context obtained from TSCreate() 98 . A - matrix 99 . B - preconditioner matrix (usually same as A) 100 . f - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 101 if A is not a function of t. 102 - ctx - [optional] user-defined context for private data for the 103 matrix evaluation routine (may be PETSC_NULL) 104 105 Calling sequence of func: 106 $ func (TS ts,double t,Mat *A,Mat *B,int *flag,void *ctx); 107 108 + t - current timestep 109 . A - matrix A, where U_t = A(t) U 110 . B - preconditioner matrix, usually the same as A 111 . flag - flag indicating information about the preconditioner matrix 112 structure (same as flag in SLESSetOperators()) 113 - ctx - [optional] user-defined context for matrix evaluation routine 114 115 Notes: 116 See SLESSetOperators() for important information about setting the flag 117 output parameter in the routine func(). Be sure to read this information! 118 119 The routine func() takes Mat * as the matrix arguments rather than Mat. 120 This allows the matrix evaluation routine to replace A and/or B with a 121 completely new new matrix structure (not just different matrix elements) 122 when appropriate, for instance, if the nonzero structure is changing 123 throughout the global iterations. 124 125 Important: 126 The user MUST call either this routine or TSSetRHSFunction(). 127 128 Level: beginner 129 130 .keywords: TS, timestep, set, right-hand-side, matrix 131 132 .seealso: TSSetRHSFunction() 133 @*/ 134 int TSSetRHSMatrix(TS ts,Mat A, Mat B,int (*f)(TS,double,Mat*,Mat*,MatStructure*,void*),void *ctx) 135 { 136 PetscFunctionBegin; 137 PetscValidHeaderSpecific(ts,TS_COOKIE); 138 if (ts->problem_type == TS_NONLINEAR) { 139 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for nonlinear problems; use TSSetRHSJacobian()"); 140 } 141 142 ts->rhsmatrix = f; 143 ts->jacP = ctx; 144 ts->A = A; 145 ts->B = B; 146 147 PetscFunctionReturn(0); 148 } 149 150 #undef __FUNC__ 151 #define __FUNC__ "TSSetRHSJacobian" 152 /*@C 153 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 154 where U_t = F(U,t), as well as the location to store the matrix. 155 156 Collective on TS 157 158 Input Parameters: 159 + ts - the TS context obtained from TSCreate() 160 . A - Jacobian matrix 161 . B - preconditioner matrix (usually same as A) 162 . f - the Jacobian evaluation routine 163 - ctx - [optional] user-defined context for private data for the 164 Jacobian evaluation routine (may be PETSC_NULL) 165 166 Calling sequence of func: 167 $ func (TS ts,double t,Vec u,Mat *A,Mat *B,int *flag,void *ctx); 168 169 + t - current timestep 170 . u - input vector 171 . A - matrix A, where U_t = A(t)u 172 . B - preconditioner matrix, usually the same as A 173 . flag - flag indicating information about the preconditioner matrix 174 structure (same as flag in SLESSetOperators()) 175 - ctx - [optional] user-defined context for matrix evaluation routine 176 177 Notes: 178 See SLESSetOperators() for important information about setting the flag 179 output parameter in the routine func(). Be sure to read this information! 180 181 The routine func() takes Mat * as the matrix arguments rather than Mat. 182 This allows the matrix evaluation routine to replace A and/or B with a 183 completely new new matrix structure (not just different matrix elements) 184 when appropriate, for instance, if the nonzero structure is changing 185 throughout the global iterations. 186 187 Level: beginner 188 189 .keywords: TS, timestep, set, right-hand-side, Jacobian 190 191 .seealso: TSDefaultComputeJacobianColor(), 192 SNESDefaultComputeJacobianColor() 193 194 @*/ 195 int TSSetRHSJacobian(TS ts,Mat A, Mat B,int (*f)(TS,double,Vec,Mat*,Mat*, 196 MatStructure*,void*),void *ctx) 197 { 198 PetscFunctionBegin; 199 PetscValidHeaderSpecific(ts,TS_COOKIE); 200 if (ts->problem_type != TS_NONLINEAR) { 201 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for linear problems; use TSSetRHSMatrix()"); 202 } 203 204 ts->rhsjacobian = f; 205 ts->jacP = ctx; 206 ts->A = A; 207 ts->B = B; 208 PetscFunctionReturn(0); 209 } 210 211 #undef __FUNC__ 212 #define __FUNC__ "TSComputeRHSBoundaryConditions" 213 /* 214 TSComputeRHSBoundaryConditions - Evaluates the boundary condition function. 215 216 Note: If the user did not provide a function but merely a matrix, 217 this routine applies the matrix. 218 */ 219 int TSComputeRHSBoundaryConditions(TS ts,double t,Vec x) 220 { 221 int ierr; 222 223 PetscFunctionBegin; 224 PetscValidHeaderSpecific(ts,TS_COOKIE); 225 PetscValidHeader(x); 226 227 if (ts->rhsbc) { 228 PetscStackPush("TS user boundary condition function"); 229 ierr = (*ts->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr); 230 PetscStackPop; 231 PetscFunctionReturn(0); 232 } 233 234 PetscFunctionReturn(0); 235 } 236 237 #undef __FUNC__ 238 #define __FUNC__ "TSSetRHSBoundaryConditions" 239 /*@C 240 TSSetRHSBoundaryConditions - Sets the routine for evaluating the function, 241 boundary conditions for the function F. 242 243 Collective on TS 244 245 Input Parameters: 246 + ts - the TS context obtained from TSCreate() 247 . f - routine for evaluating the boundary condition function 248 - ctx - [optional] user-defined context for private data for the 249 function evaluation routine (may be PETSC_NULL) 250 251 Calling sequence of func: 252 $ func (TS ts,double t,Vec F,void *ctx); 253 254 + t - current timestep 255 . F - function vector 256 - ctx - [optional] user-defined function context 257 258 Level: intermediate 259 260 .keywords: TS, timestep, set, boundary conditions, function 261 @*/ 262 int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,double,Vec,void*),void *ctx) 263 { 264 PetscFunctionBegin; 265 266 PetscValidHeaderSpecific(ts,TS_COOKIE); 267 if (ts->problem_type != TS_LINEAR) { 268 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For linear problems only"); 269 } 270 ts->rhsbc = f; 271 ts->bcP = ctx; 272 PetscFunctionReturn(0); 273 } 274 275 #undef __FUNC__ 276 #define __FUNC__ "TSView" 277 /*@ 278 TSView - Prints the TS data structure. 279 280 Collective on TS, unless Viewer is VIEWER_STDOUT_SELF 281 282 Input Parameters: 283 + ts - the TS context obtained from TSCreate() 284 - viewer - visualization context 285 286 Options Database Key: 287 . -ts_view - calls TSView() at end of TSStep() 288 289 Notes: 290 The available visualization contexts include 291 + VIEWER_STDOUT_SELF - standard output (default) 292 - VIEWER_STDOUT_WORLD - synchronized standard 293 output where only the first processor opens 294 the file. All other processors send their 295 data to the first processor to print. 296 297 The user can open an alternative visualization context with 298 ViewerASCIIOpen() - output to a specified file. 299 300 Level: beginner 301 302 .keywords: TS, timestep, view 303 304 .seealso: ViewerASCIIOpen() 305 @*/ 306 int TSView(TS ts,Viewer viewer) 307 { 308 int ierr; 309 char *method; 310 ViewerType vtype; 311 312 PetscFunctionBegin; 313 PetscValidHeaderSpecific(ts,TS_COOKIE); 314 ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 315 if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 316 ierr = ViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr); 317 ierr = TSGetType(ts,(TSType *)&method);CHKERRQ(ierr); 318 ierr = ViewerASCIIPrintf(viewer," method: %s\n",method);CHKERRQ(ierr); 319 if (ts->view) { 320 ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 321 ierr = (*ts->view)(ts,viewer);CHKERRQ(ierr); 322 ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 323 } 324 ierr = ViewerASCIIPrintf(viewer," maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr); 325 ierr = ViewerASCIIPrintf(viewer," maximum time=%g\n",ts->max_time);CHKERRQ(ierr); 326 if (ts->problem_type == TS_NONLINEAR) { 327 ierr = ViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr); 328 } 329 ierr = ViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr); 330 } else if (PetscTypeCompare(vtype,STRING_VIEWER)) { 331 ierr = TSGetType(ts,(TSType *)&method);CHKERRQ(ierr); 332 ierr = ViewerStringSPrintf(viewer," %-7.7s",method);CHKERRQ(ierr); 333 } 334 ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 335 if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);} 336 if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 337 ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 338 PetscFunctionReturn(0); 339 } 340 341 342 #undef __FUNC__ 343 #define __FUNC__ "TSSetApplicationContext" 344 /*@C 345 TSSetApplicationContext - Sets an optional user-defined context for 346 the timesteppers. 347 348 Collective on TS 349 350 Input Parameters: 351 + ts - the TS context obtained from TSCreate() 352 - usrP - optional user context 353 354 Level: intermediate 355 356 .keywords: TS, timestep, set, application, context 357 358 .seealso: TSGetApplicationContext() 359 @*/ 360 int TSSetApplicationContext(TS ts,void *usrP) 361 { 362 PetscFunctionBegin; 363 PetscValidHeaderSpecific(ts,TS_COOKIE); 364 ts->user = usrP; 365 PetscFunctionReturn(0); 366 } 367 368 #undef __FUNC__ 369 #define __FUNC__ "TSGetApplicationContext" 370 /*@C 371 TSGetApplicationContext - Gets the user-defined context for the 372 timestepper. 373 374 Not Collective 375 376 Input Parameter: 377 . ts - the TS context obtained from TSCreate() 378 379 Output Parameter: 380 . usrP - user context 381 382 Level: intermediate 383 384 .keywords: TS, timestep, get, application, context 385 386 .seealso: TSSetApplicationContext() 387 @*/ 388 int TSGetApplicationContext( TS ts, void **usrP ) 389 { 390 PetscFunctionBegin; 391 PetscValidHeaderSpecific(ts,TS_COOKIE); 392 *usrP = ts->user; 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNC__ 397 #define __FUNC__ "TSGetTimeStepNumber" 398 /*@ 399 TSGetTimeStepNumber - Gets the current number of timesteps. 400 401 Not Collective 402 403 Input Parameter: 404 . ts - the TS context obtained from TSCreate() 405 406 Output Parameter: 407 . iter - number steps so far 408 409 Level: intermediate 410 411 .keywords: TS, timestep, get, iteration, number 412 @*/ 413 int TSGetTimeStepNumber(TS ts,int* iter) 414 { 415 PetscFunctionBegin; 416 PetscValidHeaderSpecific(ts,TS_COOKIE); 417 *iter = ts->steps; 418 PetscFunctionReturn(0); 419 } 420 421 #undef __FUNC__ 422 #define __FUNC__ "TSSetInitialTimeStep" 423 /*@ 424 TSSetInitialTimeStep - Sets the initial timestep to be used, 425 as well as the initial time. 426 427 Collective on TS 428 429 Input Parameters: 430 + ts - the TS context obtained from TSCreate() 431 . initial_time - the initial time 432 - time_step - the size of the timestep 433 434 Level: intermediate 435 436 .seealso: TSSetTimeStep(), TSGetTimeStep() 437 438 .keywords: TS, set, initial, timestep 439 @*/ 440 int TSSetInitialTimeStep(TS ts,double initial_time,double time_step) 441 { 442 PetscFunctionBegin; 443 PetscValidHeaderSpecific(ts,TS_COOKIE); 444 ts->time_step = time_step; 445 ts->initial_time_step = time_step; 446 ts->ptime = initial_time; 447 PetscFunctionReturn(0); 448 } 449 450 #undef __FUNC__ 451 #define __FUNC__ "TSSetTimeStep" 452 /*@ 453 TSSetTimeStep - Allows one to reset the timestep at any time, 454 useful for simple pseudo-timestepping codes. 455 456 Collective on TS 457 458 Input Parameters: 459 + ts - the TS context obtained from TSCreate() 460 - time_step - the size of the timestep 461 462 Level: intermediate 463 464 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 465 466 .keywords: TS, set, timestep 467 @*/ 468 int TSSetTimeStep(TS ts,double time_step) 469 { 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(ts,TS_COOKIE); 472 ts->time_step = time_step; 473 PetscFunctionReturn(0); 474 } 475 476 #undef __FUNC__ 477 #define __FUNC__ "TSGetTimeStep" 478 /*@ 479 TSGetTimeStep - Gets the current timestep size. 480 481 Not Collective 482 483 Input Parameter: 484 . ts - the TS context obtained from TSCreate() 485 486 Output Parameter: 487 . dt - the current timestep size 488 489 Level: intermediate 490 491 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 492 493 .keywords: TS, get, timestep 494 @*/ 495 int TSGetTimeStep(TS ts,double* dt) 496 { 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(ts,TS_COOKIE); 499 *dt = ts->time_step; 500 PetscFunctionReturn(0); 501 } 502 503 #undef __FUNC__ 504 #define __FUNC__ "TSGetSolution" 505 /*@C 506 TSGetSolution - Returns the solution at the present timestep. It 507 is valid to call this routine inside the function that you are evaluating 508 in order to move to the new timestep. This vector not changed until 509 the solution at the next timestep has been calculated. 510 511 Not Collective, but Vec returned is parallel if TS is parallel 512 513 Input Parameter: 514 . ts - the TS context obtained from TSCreate() 515 516 Output Parameter: 517 . v - the vector containing the solution 518 519 Level: intermediate 520 521 .seealso: TSGetTimeStep() 522 523 .keywords: TS, timestep, get, solution 524 @*/ 525 int TSGetSolution(TS ts,Vec *v) 526 { 527 PetscFunctionBegin; 528 PetscValidHeaderSpecific(ts,TS_COOKIE); 529 *v = ts->vec_sol_always; 530 PetscFunctionReturn(0); 531 } 532 533 #undef __FUNC__ 534 #define __FUNC__ "TSPublish_Petsc" 535 static int TSPublish_Petsc(PetscObject object) 536 { 537 #if defined(PETSC_HAVE_AMS) 538 TS v = (TS) object; 539 int ierr; 540 #endif 541 542 PetscFunctionBegin; 543 544 #if defined(PETSC_HAVE_AMS) 545 /* if it is already published then return */ 546 if (v->amem >=0 ) PetscFunctionReturn(0); 547 548 ierr = PetscObjectPublishBaseBegin(object);CHKERRQ(ierr); 549 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ, 550 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 551 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ, 552 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 553 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1, 554 AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 555 ierr = PetscObjectPublishBaseEnd(object);CHKERRQ(ierr); 556 #endif 557 PetscFunctionReturn(0); 558 } 559 560 /* -----------------------------------------------------------*/ 561 562 #undef __FUNC__ 563 #define __FUNC__ "TSCreate" 564 /*@C 565 TSCreate - Creates a timestepper context. 566 567 Collective on MPI_Comm 568 569 Input Parameter: 570 + comm - MPI communicator 571 - type - One of TS_LINEAR,TS_NONLINEAR 572 where these types refer to problems of the forms 573 .vb 574 U_t = A U 575 U_t = A(t) U 576 U_t = F(t,U) 577 .ve 578 579 Output Parameter: 580 . outts - the new TS context 581 582 Level: beginner 583 584 .keywords: TS, timestep, create, context 585 586 .seealso: TSSetUp(), TSStep(), TSDestroy() 587 @*/ 588 int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts) 589 { 590 TS ts; 591 592 PetscFunctionBegin; 593 *outts = 0; 594 PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView); 595 PLogObjectCreate(ts); 596 ts->bops->publish = TSPublish_Petsc; 597 ts->max_steps = 5000; 598 ts->max_time = 5.0; 599 ts->time_step = .1; 600 ts->initial_time_step = ts->time_step; 601 ts->steps = 0; 602 ts->ptime = 0.0; 603 ts->data = 0; 604 ts->view = 0; 605 ts->setupcalled = 0; 606 ts->problem_type = problemtype; 607 ts->numbermonitors = 0; 608 ts->linear_its = 0; 609 ts->nonlinear_its = 0; 610 611 *outts = ts; 612 PetscFunctionReturn(0); 613 } 614 615 /* ----- Routines to initialize and destroy a timestepper ---- */ 616 617 #undef __FUNC__ 618 #define __FUNC__ "TSSetUp" 619 /*@ 620 TSSetUp - Sets up the internal data structures for the later use 621 of a timestepper. 622 623 Collective on TS 624 625 Input Parameter: 626 . ts - the TS context obtained from TSCreate() 627 628 Notes: 629 For basic use of the TS solvers the user need not explicitly call 630 TSSetUp(), since these actions will automatically occur during 631 the call to TSStep(). However, if one wishes to control this 632 phase separately, TSSetUp() should be called after TSCreate() 633 and optional routines of the form TSSetXXX(), but before TSStep(). 634 635 Level: advanced 636 637 .keywords: TS, timestep, setup 638 639 .seealso: TSCreate(), TSStep(), TSDestroy() 640 @*/ 641 int TSSetUp(TS ts) 642 { 643 int ierr; 644 645 PetscFunctionBegin; 646 PetscValidHeaderSpecific(ts,TS_COOKIE); 647 if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call TSSetSolution() first"); 648 if (!ts->type_name) { 649 ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr); 650 } 651 ierr = (*ts->setup)(ts);CHKERRQ(ierr); 652 ts->setupcalled = 1; 653 PetscFunctionReturn(0); 654 } 655 656 #undef __FUNC__ 657 #define __FUNC__ "TSDestroy" 658 /*@C 659 TSDestroy - Destroys the timestepper context that was created 660 with TSCreate(). 661 662 Collective on TS 663 664 Input Parameter: 665 . ts - the TS context obtained from TSCreate() 666 667 Level: beginner 668 669 .keywords: TS, timestepper, destroy 670 671 .seealso: TSCreate(), TSSetUp(), TSSolve() 672 @*/ 673 int TSDestroy(TS ts) 674 { 675 int ierr; 676 677 PetscFunctionBegin; 678 PetscValidHeaderSpecific(ts,TS_COOKIE); 679 if (--ts->refct > 0) PetscFunctionReturn(0); 680 681 /* if memory was published with AMS then destroy it */ 682 ierr = PetscAMSDestroy(ts);CHKERRQ(ierr); 683 684 if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);} 685 if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 686 ierr = (*(ts)->destroy)(ts);CHKERRQ(ierr); 687 PLogObjectDestroy((PetscObject)ts); 688 PetscHeaderDestroy((PetscObject)ts); 689 PetscFunctionReturn(0); 690 } 691 692 #undef __FUNC__ 693 #define __FUNC__ "TSGetSNES" 694 /*@C 695 TSGetSNES - Returns the SNES (nonlinear solver) associated with 696 a TS (timestepper) context. Valid only for nonlinear problems. 697 698 Not Collective, but SNES is parallel if TS is parallel 699 700 Input Parameter: 701 . ts - the TS context obtained from TSCreate() 702 703 Output Parameter: 704 . snes - the nonlinear solver context 705 706 Notes: 707 The user can then directly manipulate the SNES context to set various 708 options, etc. Likewise, the user can then extract and manipulate the 709 SLES, KSP, and PC contexts as well. 710 711 TSGetSNES() does not work for integrators that do not use SNES; in 712 this case TSGetSNES() returns PETSC_NULL in snes. 713 714 Level: beginner 715 716 .keywords: timestep, get, SNES 717 @*/ 718 int TSGetSNES(TS ts,SNES *snes) 719 { 720 PetscFunctionBegin; 721 PetscValidHeaderSpecific(ts,TS_COOKIE); 722 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Nonlinear only; use TSGetSLES()"); 723 *snes = ts->snes; 724 PetscFunctionReturn(0); 725 } 726 727 #undef __FUNC__ 728 #define __FUNC__ "TSGetSLES" 729 /*@C 730 TSGetSLES - Returns the SLES (linear solver) associated with 731 a TS (timestepper) context. 732 733 Not Collective, but SLES is parallel if TS is parallel 734 735 Input Parameter: 736 . ts - the TS context obtained from TSCreate() 737 738 Output Parameter: 739 . sles - the nonlinear solver context 740 741 Notes: 742 The user can then directly manipulate the SLES context to set various 743 options, etc. Likewise, the user can then extract and manipulate the 744 KSP and PC contexts as well. 745 746 TSGetSLES() does not work for integrators that do not use SLES; 747 in this case TSGetSLES() returns PETSC_NULL in sles. 748 749 Level: beginner 750 751 .keywords: timestep, get, SLES 752 @*/ 753 int TSGetSLES(TS ts,SLES *sles) 754 { 755 PetscFunctionBegin; 756 PetscValidHeaderSpecific(ts,TS_COOKIE); 757 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Linear only; use TSGetSNES()"); 758 *sles = ts->sles; 759 PetscFunctionReturn(0); 760 } 761 762 /* ----------- Routines to set solver parameters ---------- */ 763 764 #undef __FUNC__ 765 #define __FUNC__ "TSSetDuration" 766 /*@ 767 TSSetDuration - Sets the maximum number of timesteps to use and 768 maximum time for iteration. 769 770 Collective on TS 771 772 Input Parameters: 773 + ts - the TS context obtained from TSCreate() 774 . maxsteps - maximum number of iterations to use 775 - maxtime - final time to iterate to 776 777 Options Database Keys: 778 . -ts_max_steps <maxsteps> - Sets maxsteps 779 . -ts_max_time <maxtime> - Sets maxtime 780 781 Notes: 782 The default maximum number of iterations is 5000. Default time is 5.0 783 784 Level: intermediate 785 786 .keywords: TS, timestep, set, maximum, iterations 787 @*/ 788 int TSSetDuration(TS ts,int maxsteps,double maxtime) 789 { 790 PetscFunctionBegin; 791 PetscValidHeaderSpecific(ts,TS_COOKIE); 792 ts->max_steps = maxsteps; 793 ts->max_time = maxtime; 794 PetscFunctionReturn(0); 795 } 796 797 #undef __FUNC__ 798 #define __FUNC__ "TSSetSolution" 799 /*@ 800 TSSetSolution - Sets the initial solution vector 801 for use by the TS routines. 802 803 Collective on TS and Vec 804 805 Input Parameters: 806 + ts - the TS context obtained from TSCreate() 807 - x - the solution vector 808 809 Level: beginner 810 811 .keywords: TS, timestep, set, solution, initial conditions 812 @*/ 813 int TSSetSolution(TS ts,Vec x) 814 { 815 PetscFunctionBegin; 816 PetscValidHeaderSpecific(ts,TS_COOKIE); 817 ts->vec_sol = ts->vec_sol_always = x; 818 PetscFunctionReturn(0); 819 } 820 821 /* ------------ Routines to set performance monitoring options ----------- */ 822 823 #undef __FUNC__ 824 #define __FUNC__ "TSSetMonitor" 825 /*@C 826 TSSetMonitor - Sets an ADDITIONAL function that is to be used at every 827 timestep to display the iteration's progress. 828 829 Collective on TS 830 831 Input Parameters: 832 + ts - the TS context obtained from TSCreate() 833 . func - monitoring routine 834 - mctx - [optional] user-defined context for private data for the 835 monitor routine (may be PETSC_NULL) 836 837 Calling sequence of func: 838 $ int func(TS ts,int steps,double time,Vec x,void *mctx) 839 840 + ts - the TS context 841 . steps - iteration number 842 . time - current timestep 843 . x - current iterate 844 - mctx - [optional] monitoring context 845 846 Notes: 847 This routine adds an additional monitor to the list of monitors that 848 already has been loaded. 849 850 Level: intermediate 851 852 .keywords: TS, timestep, set, monitor 853 854 .seealso: TSDefaultMonitor(), TSClearMonitor() 855 @*/ 856 int TSSetMonitor(TS ts, int (*monitor)(TS,int,double,Vec,void*), void *mctx ) 857 { 858 PetscFunctionBegin; 859 PetscValidHeaderSpecific(ts,TS_COOKIE); 860 if (ts->numbermonitors >= MAXTSMONITORS) { 861 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 862 } 863 ts->monitor[ts->numbermonitors] = monitor; 864 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 865 PetscFunctionReturn(0); 866 } 867 868 #undef __FUNC__ 869 #define __FUNC__ "TSClearMonitor" 870 /*@C 871 TSClearMonitor - Clears all the monitors that have been set on a time-step object. 872 873 Collective on TS 874 875 Input Parameters: 876 . ts - the TS context obtained from TSCreate() 877 878 Notes: 879 There is no way to remove a single, specific monitor. 880 881 Level: intermediate 882 883 .keywords: TS, timestep, set, monitor 884 885 .seealso: TSDefaultMonitor(), TSSetMonitor() 886 @*/ 887 int TSClearMonitor(TS ts) 888 { 889 PetscFunctionBegin; 890 PetscValidHeaderSpecific(ts,TS_COOKIE); 891 ts->numbermonitors = 0; 892 PetscFunctionReturn(0); 893 } 894 895 #undef __FUNC__ 896 #define __FUNC__ "TSDefaultMonitor" 897 int TSDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 898 { 899 int ierr; 900 901 PetscFunctionBegin; 902 ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,time);CHKERRQ(ierr); 903 PetscFunctionReturn(0); 904 } 905 906 #undef __FUNC__ 907 #define __FUNC__ "TSStep" 908 /*@ 909 TSStep - Steps the requested number of timesteps. 910 911 Collective on TS 912 913 Input Parameter: 914 . ts - the TS context obtained from TSCreate() 915 916 Output Parameters: 917 + steps - number of iterations until termination 918 - time - time until termination 919 920 Level: beginner 921 922 .keywords: TS, timestep, solve 923 924 .seealso: TSCreate(), TSSetUp(), TSDestroy() 925 @*/ 926 int TSStep(TS ts,int *steps,double *time) 927 { 928 int ierr,flg; 929 930 PetscFunctionBegin; 931 PetscValidHeaderSpecific(ts,TS_COOKIE); 932 if (!ts->setupcalled) {ierr = TSSetUp(ts);CHKERRQ(ierr);} 933 PLogEventBegin(TS_Step,ts,0,0,0); 934 ierr = (*(ts)->step)(ts,steps,time);CHKERRQ(ierr); 935 PLogEventEnd(TS_Step,ts,0,0,0); 936 ierr = OptionsHasName(PETSC_NULL,"-ts_view",&flg);CHKERRQ(ierr); 937 if (flg) { 938 ierr = TSView(ts,VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 939 } 940 PetscFunctionReturn(0); 941 } 942 943 #undef __FUNC__ 944 #define __FUNC__ "TSMonitor" 945 /* 946 Runs the user provided monitor routines, if they exists. 947 */ 948 int TSMonitor(TS ts,int step,double time,Vec x) 949 { 950 int i,ierr,n = ts->numbermonitors; 951 952 PetscFunctionBegin; 953 for ( i=0; i<n; i++ ) { 954 ierr = (*ts->monitor[i])(ts,step,time,x,ts->monitorcontext[i]);CHKERRQ(ierr); 955 } 956 PetscFunctionReturn(0); 957 } 958 959 /* ------------------------------------------------------------------------*/ 960 961 #undef __FUNC__ 962 #define __FUNC__ "TSLGMonitorCreate" 963 /*@C 964 TSLGMonitorCreate - Creates a line graph context for use with 965 TS to monitor convergence of preconditioned residual norms. 966 967 Collective on TS 968 969 Input Parameters: 970 + host - the X display to open, or null for the local machine 971 . label - the title to put in the title bar 972 . x, y - the screen coordinates of the upper left coordinate of 973 the window 974 - m, n - the screen width and height in pixels 975 976 Output Parameter: 977 . draw - the drawing context 978 979 Options Database Key: 980 . -ts_xmonitor - automatically sets line graph monitor 981 982 Notes: 983 Use TSLGMonitorDestroy() to destroy this line graph, not DrawLGDestroy(). 984 985 Level: intermediate 986 987 .keywords: TS, monitor, line graph, residual, create 988 989 .seealso: TSLGMonitorDestroy(), TSSetMonitor() 990 @*/ 991 int TSLGMonitorCreate(char *host,char *label,int x,int y,int m, 992 int n, DrawLG *draw) 993 { 994 Draw win; 995 int ierr; 996 997 PetscFunctionBegin; 998 ierr = DrawOpenX(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 999 ierr = DrawLGCreate(win,1,draw);CHKERRQ(ierr); 1000 ierr = DrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1001 1002 PLogObjectParent(*draw,win); 1003 PetscFunctionReturn(0); 1004 } 1005 1006 #undef __FUNC__ 1007 #define __FUNC__ "TSLGMonitor" 1008 int TSLGMonitor(TS ts,int n,double time,Vec v,void *monctx) 1009 { 1010 DrawLG lg = (DrawLG) monctx; 1011 double x,y = time; 1012 int ierr; 1013 1014 PetscFunctionBegin; 1015 if (!n) {ierr = DrawLGReset(lg);CHKERRQ(ierr);} 1016 x = (double) n; 1017 ierr = DrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1018 if (n < 20 || (n % 5)) { 1019 ierr = DrawLGDraw(lg);CHKERRQ(ierr); 1020 } 1021 PetscFunctionReturn(0); 1022 } 1023 1024 #undef __FUNC__ 1025 #define __FUNC__ "TSLGMonitorDestroy" 1026 /*@C 1027 TSLGMonitorDestroy - Destroys a line graph context that was created 1028 with TSLGMonitorCreate(). 1029 1030 Collective on DrawLG 1031 1032 Input Parameter: 1033 . draw - the drawing context 1034 1035 Level: intermediate 1036 1037 .keywords: TS, monitor, line graph, destroy 1038 1039 .seealso: TSLGMonitorCreate(), TSSetMonitor(), TSLGMonitor(); 1040 @*/ 1041 int TSLGMonitorDestroy(DrawLG drawlg) 1042 { 1043 Draw draw; 1044 int ierr; 1045 1046 PetscFunctionBegin; 1047 ierr = DrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1048 ierr = DrawDestroy(draw);CHKERRQ(ierr); 1049 ierr = DrawLGDestroy(drawlg);CHKERRQ(ierr); 1050 PetscFunctionReturn(0); 1051 } 1052 1053 #undef __FUNC__ 1054 #define __FUNC__ "TSGetTime" 1055 /*@ 1056 TSGetTime - Gets the current time. 1057 1058 Not Collective 1059 1060 Input Parameter: 1061 . ts - the TS context obtained from TSCreate() 1062 1063 Output Parameter: 1064 . t - the current time 1065 1066 Contributed by: Matthew Knepley 1067 1068 Level: beginner 1069 1070 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1071 1072 .keywords: TS, get, time 1073 @*/ 1074 int TSGetTime(TS ts, double* t) 1075 { 1076 PetscFunctionBegin; 1077 PetscValidHeaderSpecific(ts, TS_COOKIE); 1078 *t = ts->ptime; 1079 PetscFunctionReturn(0); 1080 } 1081 1082 #undef __FUNC__ 1083 #define __FUNC__ "TSGetProblemType" 1084 /*@C 1085 TSGetProblemType - Returns the problem type of a TS (timestepper) context. 1086 1087 Not Collective 1088 1089 Input Parameter: 1090 . ts - The TS context obtained from TSCreate() 1091 1092 Output Parameter: 1093 . type - The problem type, TS_LINEAR or TS_NONLINEAR 1094 1095 Level: intermediate 1096 1097 Contributed by: Matthew Knepley 1098 1099 .keywords: ts, get, type 1100 1101 @*/ 1102 int TSGetProblemType(TS ts, TSProblemType *type) 1103 { 1104 PetscFunctionBegin; 1105 PetscValidHeaderSpecific(ts, TS_COOKIE); 1106 *type = ts->problem_type; 1107 PetscFunctionReturn(0); 1108 } 1109 1110 #undef __FUNC__ 1111 #define __FUNC__ "TSSetOptionsPrefix" 1112 /*@C 1113 TSSetOptionsPrefix - Sets the prefix used for searching for all 1114 TS options in the database. 1115 1116 Collective on TS 1117 1118 Input Parameter: 1119 + ts - The TS context 1120 - prefix - The prefix to prepend to all option names 1121 1122 Notes: 1123 A hyphen (-) must NOT be given at the beginning of the prefix name. 1124 The first character of all runtime options is AUTOMATICALLY the 1125 hyphen. 1126 1127 Contributed by: Matthew Knepley 1128 1129 Level: advanced 1130 1131 .keywords: TS, set, options, prefix, database 1132 1133 .seealso: TSSetFromOptions() 1134 1135 @*/ 1136 int TSSetOptionsPrefix(TS ts, char *prefix) 1137 { 1138 int ierr; 1139 1140 PetscFunctionBegin; 1141 PetscValidHeaderSpecific(ts, TS_COOKIE); 1142 ierr = PetscObjectSetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr); 1143 switch(ts->problem_type) { 1144 case TS_NONLINEAR: 1145 ierr = SNESSetOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr); 1146 break; 1147 case TS_LINEAR: 1148 ierr = SLESSetOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr); 1149 break; 1150 } 1151 PetscFunctionReturn(0); 1152 } 1153 1154 1155 #undef __FUNC__ 1156 #define __FUNC__ "TSAppendOptionsPrefix" 1157 /*@C 1158 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1159 TS options in the database. 1160 1161 Collective on TS 1162 1163 Input Parameter: 1164 + ts - The TS context 1165 - prefix - The prefix to prepend to all option names 1166 1167 Notes: 1168 A hyphen (-) must NOT be given at the beginning of the prefix name. 1169 The first character of all runtime options is AUTOMATICALLY the 1170 hyphen. 1171 1172 Contributed by: Matthew Knepley 1173 1174 Level: advanced 1175 1176 .keywords: TS, append, options, prefix, database 1177 1178 .seealso: TSGetOptionsPrefix() 1179 1180 @*/ 1181 int TSAppendOptionsPrefix(TS ts, char *prefix) 1182 { 1183 int ierr; 1184 1185 PetscFunctionBegin; 1186 PetscValidHeaderSpecific(ts, TS_COOKIE); 1187 ierr = PetscObjectAppendOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr); 1188 switch(ts->problem_type) { 1189 case TS_NONLINEAR: 1190 ierr = SNESAppendOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr); 1191 break; 1192 case TS_LINEAR: 1193 ierr = SLESAppendOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr); 1194 break; 1195 } 1196 PetscFunctionReturn(0); 1197 } 1198 1199 #undef __FUNC__ 1200 #define __FUNC__ "TSGetOptionsPrefix" 1201 /*@C 1202 TSGetOptionsPrefix - Sets the prefix used for searching for all 1203 TS options in the database. 1204 1205 Not Collective 1206 1207 Input Parameter: 1208 . ts - The TS context 1209 1210 Output Parameter: 1211 . prefix - A pointer to the prefix string used 1212 1213 Contributed by: Matthew Knepley 1214 1215 Notes: On the fortran side, the user should pass in a string 'prifix' of 1216 sufficient length to hold the prefix. 1217 1218 Level: intermediate 1219 1220 .keywords: TS, get, options, prefix, database 1221 1222 .seealso: TSAppendOptionsPrefix() 1223 @*/ 1224 int TSGetOptionsPrefix(TS ts, char **prefix) 1225 { 1226 int ierr; 1227 1228 PetscFunctionBegin; 1229 PetscValidHeaderSpecific(ts, TS_COOKIE); 1230 ierr = PetscObjectGetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr); 1231 PetscFunctionReturn(0); 1232 } 1233 1234 #undef __FUNC__ 1235 #define __FUNC__ "TSGetRHSMatrix" 1236 /*@C 1237 TSGetRHSMatrix - Returns the matrix A at the present timestep. 1238 1239 Not Collective, but parallel objects are returned if TS is parallel 1240 1241 Input Parameter: 1242 . ts - The TS context obtained from TSCreate() 1243 1244 Output Parameters: 1245 + A - The matrix A, where U_t = A(t) U 1246 . M - The preconditioner matrix, usually the same as A 1247 - ctx - User-defined context for matrix evaluation routine 1248 1249 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1250 1251 Contributed by: Matthew Knepley 1252 1253 Level: intermediate 1254 1255 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 1256 1257 .keywords: TS, timestep, get, matrix 1258 1259 @*/ 1260 int TSGetRHSMatrix(TS ts, Mat *A, Mat *M, void **ctx) 1261 { 1262 PetscFunctionBegin; 1263 PetscValidHeaderSpecific(ts, TS_COOKIE); 1264 if (A) *A = ts->A; 1265 if (M) *M = ts->B; 1266 if (ctx) *ctx = ts->jacP; 1267 PetscFunctionReturn(0); 1268 } 1269 1270 #undef __FUNC__ 1271 #define __FUNC__ "TSGetRHSJacobian" 1272 /*@C 1273 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1274 1275 Not Collective, but parallel objects are returned if TS is parallel 1276 1277 Input Parameter: 1278 . ts - The TS context obtained from TSCreate() 1279 1280 Output Parameters: 1281 + J - The Jacobian J of F, where U_t = F(U,t) 1282 . M - The preconditioner matrix, usually the same as J 1283 - ctx - User-defined context for Jacobian evaluation routine 1284 1285 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1286 1287 Contributed by: Matthew Knepley 1288 1289 Level: intermediate 1290 1291 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber() 1292 1293 .keywords: TS, timestep, get, matrix, Jacobian 1294 @*/ 1295 int TSGetRHSJacobian(TS ts, Mat *J, Mat *M, void **ctx) 1296 { 1297 int ierr; 1298 1299 PetscFunctionBegin; 1300 ierr = TSGetRHSMatrix(ts, J, M, ctx);CHKERRQ(ierr); 1301 PetscFunctionReturn(0); 1302 } 1303 1304 /*MC 1305 TSRegister - Adds a method to the timestepping solver package. 1306 1307 Synopsis: 1308 1309 TSRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(TS)) 1310 1311 Not collective 1312 1313 Input Parameters: 1314 + name_solver - name of a new user-defined solver 1315 . path - path (either absolute or relative) the library containing this solver 1316 . name_create - name of routine to create method context 1317 - routine_create - routine to create method context 1318 1319 Notes: 1320 TSRegister() may be called multiple times to add several user-defined solvers. 1321 1322 If dynamic libraries are used, then the fourth input argument (routine_create) 1323 is ignored. 1324 1325 Sample usage: 1326 .vb 1327 TSRegister("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 1328 "MySolverCreate",MySolverCreate); 1329 .ve 1330 1331 Then, your solver can be chosen with the procedural interface via 1332 $ TSSetType(ts,"my_solver") 1333 or at runtime via the option 1334 $ -ts_type my_solver 1335 1336 Level: advanced 1337 1338 $PETSC_ARCH, $PETSC_DIR, $PETSC_LDIR, and $BOPT occuring in pathname will be replaced with appropriate values. 1339 1340 .keywords: TS, register 1341 1342 .seealso: TSRegisterAll(), TSRegisterDestroy() 1343 M*/ 1344 1345 #undef __FUNC__ 1346 #define __FUNC__ "TSRegister_Private" 1347 int TSRegister_Private(char *sname,char *path,char *name,int (*function)(TS)) 1348 { 1349 char fullname[256]; 1350 int ierr; 1351 1352 PetscFunctionBegin; 1353 ierr = PetscStrcpy(fullname,path);CHKERRQ(ierr); 1354 PetscStrcat(fullname,":"); PetscStrcat(fullname,name); 1355 FListAdd_Private(&TSList,sname,fullname, (int (*)(void*))function); 1356 PetscFunctionReturn(0); 1357 } 1358