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