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