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