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->ops->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->ops->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->ops->rhsfunction) { 88 PetscStackPush("TS user right-hand-side function"); 89 ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 90 PetscStackPop; 91 } else { 92 if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 93 MatStructure flg; 94 PetscStackPush("TS user right-hand-side matrix function"); 95 ierr = (*ts->ops->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->ops->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->ops->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->ops->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->ops->rhsbc) { 300 PetscStackPush("TS user boundary condition function"); 301 ierr = (*ts->ops->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->ops->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->ops->view) { 401 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 402 ierr = (*ts->ops->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 . ts - 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 *ts) 670 { 671 TS t; 672 int ierr; 673 674 PetscFunctionBegin; 675 *ts = PETSC_NULL; 676 PetscValidPointer(ts); 677 PetscHeaderCreate(t, _p_TS, struct _TSOps, TS_COOKIE, -1, "TS", comm, TSDestroy, TSView); 678 PetscLogObjectCreate(t); 679 PetscLogObjectMemory(t, sizeof(struct _p_TS)); 680 ierr = PetscMemzero(t->ops, sizeof(struct _TSOps)); CHKERRQ(ierr); 681 t->bops->publish = TSPublish_Petsc; 682 t->type_name = PETSC_NULL; 683 684 t->problem_type = problemtype; 685 t->vec_sol = PETSC_NULL; 686 t->vec_sol_always = PETSC_NULL; 687 t->numbermonitors = 0; 688 t->isGTS = PETSC_FALSE; 689 t->explicit = PETSC_NULL; 690 t->Iindex = PETSC_NULL; 691 t->sles = PETSC_NULL; 692 t->A = PETSC_NULL; 693 t->B = PETSC_NULL; 694 t->snes = PETSC_NULL; 695 t->funP = PETSC_NULL; 696 t->jacP = PETSC_NULL; 697 t->setupcalled = 0; 698 t->data = PETSC_NULL; 699 t->user = PETSC_NULL; 700 t->max_steps = 5000; 701 t->max_time = 5.0; 702 t->time_step = .1; 703 t->time_step_old = t->time_step; 704 t->initial_time_step = t->time_step; 705 t->steps = 0; 706 t->ptime = 0.0; 707 t->linear_its = 0; 708 t->nonlinear_its = 0; 709 t->work = PETSC_NULL; 710 t->nwork = 0; 711 t->ops->applymatrixbc = TSDefaultSystemMatrixBC; 712 t->ops->applyrhsbc = TSDefaultRhsBC; 713 t->ops->applysolbc = TSDefaultSolutionBC; 714 t->ops->prestep = TSDefaultPreStep; 715 t->ops->update = TSDefaultUpdate; 716 t->ops->poststep = TSDefaultPostStep; 717 t->ops->reform = PETSC_NULL; 718 t->ops->reallocate = PETSC_NULL; 719 t->ops->serialize = PETSC_NULL; 720 721 *ts = t; 722 PetscFunctionReturn(0); 723 } 724 725 /* ----- Routines to initialize and destroy a timestepper ---- */ 726 727 #undef __FUNCT__ 728 #define __FUNCT__ "TSSetUp" 729 /*@ 730 TSSetUp - Sets up the internal data structures for the later use 731 of a timestepper. 732 733 Collective on TS 734 735 Input Parameter: 736 . ts - the TS context obtained from TSCreate() 737 738 Notes: 739 For basic use of the TS solvers the user need not explicitly call 740 TSSetUp(), since these actions will automatically occur during 741 the call to TSStep(). However, if one wishes to control this 742 phase separately, TSSetUp() should be called after TSCreate() 743 and optional routines of the form TSSetXXX(), but before TSStep(). 744 745 Level: advanced 746 747 .keywords: TS, timestep, setup 748 749 .seealso: TSCreate(), TSStep(), TSDestroy() 750 @*/ 751 int TSSetUp(TS ts) 752 { 753 int ierr; 754 755 PetscFunctionBegin; 756 PetscValidHeaderSpecific(ts,TS_COOKIE); 757 if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 758 if (!ts->type_name) { 759 ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr); 760 } 761 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 762 ts->setupcalled = 1; 763 PetscFunctionReturn(0); 764 } 765 766 #undef __FUNCT__ 767 #define __FUNCT__ "TSDestroy" 768 /*@C 769 TSDestroy - Destroys the timestepper context that was created 770 with TSCreate(). 771 772 Collective on TS 773 774 Input Parameter: 775 . ts - the TS context obtained from TSCreate() 776 777 Level: beginner 778 779 .keywords: TS, timestepper, destroy 780 781 .seealso: TSCreate(), TSSetUp(), TSSolve() 782 @*/ 783 int TSDestroy(TS ts) 784 { 785 int ierr,i; 786 787 PetscFunctionBegin; 788 PetscValidHeaderSpecific(ts,TS_COOKIE); 789 if (--ts->refct > 0) PetscFunctionReturn(0); 790 791 /* if memory was published with AMS then destroy it */ 792 ierr = PetscObjectDepublish(ts);CHKERRQ(ierr); 793 794 if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);} 795 if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 796 ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr); 797 for (i=0; i<ts->numbermonitors; i++) { 798 if (ts->mdestroy[i]) { 799 ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr); 800 } 801 } 802 PetscLogObjectDestroy((PetscObject)ts); 803 PetscHeaderDestroy((PetscObject)ts); 804 PetscFunctionReturn(0); 805 } 806 807 #undef __FUNCT__ 808 #define __FUNCT__ "TSGetSNES" 809 /*@C 810 TSGetSNES - Returns the SNES (nonlinear solver) associated with 811 a TS (timestepper) context. Valid only for nonlinear problems. 812 813 Not Collective, but SNES is parallel if TS is parallel 814 815 Input Parameter: 816 . ts - the TS context obtained from TSCreate() 817 818 Output Parameter: 819 . snes - the nonlinear solver context 820 821 Notes: 822 The user can then directly manipulate the SNES context to set various 823 options, etc. Likewise, the user can then extract and manipulate the 824 SLES, KSP, and PC contexts as well. 825 826 TSGetSNES() does not work for integrators that do not use SNES; in 827 this case TSGetSNES() returns PETSC_NULL in snes. 828 829 Level: beginner 830 831 .keywords: timestep, get, SNES 832 @*/ 833 int TSGetSNES(TS ts,SNES *snes) 834 { 835 PetscFunctionBegin; 836 PetscValidHeaderSpecific(ts,TS_COOKIE); 837 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetSLES()"); 838 *snes = ts->snes; 839 PetscFunctionReturn(0); 840 } 841 842 #undef __FUNCT__ 843 #define __FUNCT__ "TSGetSLES" 844 /*@C 845 TSGetSLES - Returns the SLES (linear solver) associated with 846 a TS (timestepper) context. 847 848 Not Collective, but SLES is parallel if TS is parallel 849 850 Input Parameter: 851 . ts - the TS context obtained from TSCreate() 852 853 Output Parameter: 854 . sles - the nonlinear solver context 855 856 Notes: 857 The user can then directly manipulate the SLES context to set various 858 options, etc. Likewise, the user can then extract and manipulate the 859 KSP and PC contexts as well. 860 861 TSGetSLES() does not work for integrators that do not use SLES; 862 in this case TSGetSLES() returns PETSC_NULL in sles. 863 864 Level: beginner 865 866 .keywords: timestep, get, SLES 867 @*/ 868 int TSGetSLES(TS ts,SLES *sles) 869 { 870 PetscFunctionBegin; 871 PetscValidHeaderSpecific(ts,TS_COOKIE); 872 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 873 *sles = ts->sles; 874 PetscFunctionReturn(0); 875 } 876 877 /* ----------- Routines to set solver parameters ---------- */ 878 879 #undef __FUNCT__ 880 #define __FUNCT__ "TSSetDuration" 881 /*@ 882 TSSetDuration - Sets the maximum number of timesteps to use and 883 maximum time for iteration. 884 885 Collective on TS 886 887 Input Parameters: 888 + ts - the TS context obtained from TSCreate() 889 . maxsteps - maximum number of iterations to use 890 - maxtime - final time to iterate to 891 892 Options Database Keys: 893 . -ts_max_steps <maxsteps> - Sets maxsteps 894 . -ts_max_time <maxtime> - Sets maxtime 895 896 Notes: 897 The default maximum number of iterations is 5000. Default time is 5.0 898 899 Level: intermediate 900 901 .keywords: TS, timestep, set, maximum, iterations 902 @*/ 903 int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime) 904 { 905 PetscFunctionBegin; 906 PetscValidHeaderSpecific(ts,TS_COOKIE); 907 ts->max_steps = maxsteps; 908 ts->max_time = maxtime; 909 PetscFunctionReturn(0); 910 } 911 912 #undef __FUNCT__ 913 #define __FUNCT__ "TSSetSolution" 914 /*@ 915 TSSetSolution - Sets the initial solution vector 916 for use by the TS routines. 917 918 Collective on TS and Vec 919 920 Input Parameters: 921 + ts - the TS context obtained from TSCreate() 922 - x - the solution vector 923 924 Level: beginner 925 926 .keywords: TS, timestep, set, solution, initial conditions 927 @*/ 928 int TSSetSolution(TS ts,Vec x) 929 { 930 PetscFunctionBegin; 931 PetscValidHeaderSpecific(ts,TS_COOKIE); 932 ts->vec_sol = ts->vec_sol_always = x; 933 PetscFunctionReturn(0); 934 } 935 936 #undef __FUNC__ 937 #define __FUNC__ "TSSetRhsBC" 938 /*@ 939 TSSetRhsBC - Sets the function which applies boundary conditions 940 to the Rhs of each system. 941 942 Collective on TS 943 944 Input Parameters: 945 + ts - The TS context obtained from TSCreate() 946 - func - The function 947 948 Calling sequence of func: 949 . func (TS ts, Vec rhs, void *ctx); 950 951 + rhs - The current rhs vector 952 - ctx - The user-context 953 954 Level: intermediate 955 956 .keywords: TS, Rhs, boundary conditions 957 @*/ 958 int TSSetRhsBC(TS ts, int (*func)(TS, Vec, void *)) 959 { 960 PetscFunctionBegin; 961 PetscValidHeaderSpecific(ts, TS_COOKIE); 962 ts->ops->applyrhsbc = func; 963 PetscFunctionReturn(0); 964 } 965 966 #undef __FUNC__ 967 #define __FUNC__ "TSDefaultRhsBC" 968 /*@ 969 TSDefaultRhsBC - The default boundary condition function which does nothing. 970 971 Collective on TS 972 973 Input Parameters: 974 + ts - The TS context obtained from TSCreate() 975 . rhs - The Rhs 976 - ctx - The user-context 977 978 Level: developer 979 980 .keywords: TS, Rhs, boundary conditions 981 @*/ 982 int TSDefaultRhsBC(TS ts, Vec rhs, void *ctx) 983 { 984 PetscFunctionBegin; 985 PetscFunctionReturn(0); 986 } 987 988 #undef __FUNC__ 989 #define __FUNC__ "TSSetSystemMatrixBC" 990 /*@ 991 TSSetSystemMatrixBC - Sets the function which applies boundary conditions 992 to the system matrix and preconditioner of each system. 993 994 Collective on TS 995 996 Input Parameters: 997 + ts - The TS context obtained from TSCreate() 998 - func - The function 999 1000 Calling sequence of func: 1001 . func (TS ts, Mat A, Mat B, void *ctx); 1002 1003 + A - The current system matrix 1004 . B - The current preconditioner 1005 - ctx - The user-context 1006 1007 Level: intermediate 1008 1009 .keywords: TS, System matrix, boundary conditions 1010 @*/ 1011 int TSSetSystemMatrixBC(TS ts, int (*func)(TS, Mat, Mat, void *)) 1012 { 1013 PetscFunctionBegin; 1014 PetscValidHeaderSpecific(ts, TS_COOKIE); 1015 ts->ops->applymatrixbc = func; 1016 PetscFunctionReturn(0); 1017 } 1018 1019 #undef __FUNC__ 1020 #define __FUNC__ "TSDefaultSystemMatrixBC" 1021 /*@ 1022 TSDefaultSystemMatrixBC - The default boundary condition function which 1023 does nothing. 1024 1025 Collective on TS 1026 1027 Input Parameters: 1028 + ts - The TS context obtained from TSCreate() 1029 . A - The system matrix 1030 . B - The preconditioner 1031 - ctx - The user-context 1032 1033 Level: developer 1034 1035 .keywords: TS, System matrix, boundary conditions 1036 @*/ 1037 int TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx) 1038 { 1039 PetscFunctionBegin; 1040 PetscFunctionReturn(0); 1041 } 1042 1043 #undef __FUNC__ 1044 #define __FUNC__ "TSSetSolutionBC" 1045 /*@ 1046 TSSetSolutionBC - Sets the function which applies boundary conditions 1047 to the solution of each system. This is necessary in nonlinear systems 1048 which time dependent boundary conditions. 1049 1050 Collective on TS 1051 1052 Input Parameters: 1053 + ts - The TS context obtained from TSCreate() 1054 - func - The function 1055 1056 Calling sequence of func: 1057 . func (TS ts, Vec rsol, void *ctx); 1058 1059 + sol - The current solution vector 1060 - ctx - The user-context 1061 1062 Level: intermediate 1063 1064 .keywords: TS, solution, boundary conditions 1065 @*/ 1066 int TSSetSolutionBC(TS ts, int (*func)(TS, Vec, void *)) 1067 { 1068 PetscFunctionBegin; 1069 PetscValidHeaderSpecific(ts, TS_COOKIE); 1070 ts->ops->applysolbc = func; 1071 PetscFunctionReturn(0); 1072 } 1073 1074 #undef __FUNC__ 1075 #define __FUNC__ "TSDefaultSolutionBC" 1076 /*@ 1077 TSDefaultSolutionBC - The default boundary condition function which 1078 does nothing. 1079 1080 Collective on TS 1081 1082 Input Parameters: 1083 + ts - The TS context obtained from TSCreate() 1084 . sol - The solution 1085 - ctx - The user-context 1086 1087 Level: developer 1088 1089 .keywords: TS, solution, boundary conditions 1090 @*/ 1091 int TSDefaultSolutionBC(TS ts, Vec sol, void *ctx) 1092 { 1093 PetscFunctionBegin; 1094 PetscFunctionReturn(0); 1095 } 1096 1097 #undef __FUNC__ 1098 #define __FUNC__ "TSSetPreStep" 1099 /*@ 1100 TSSetPreStep - Sets the general-purpose function 1101 called once at the beginning of time stepping. 1102 1103 Collective on TS 1104 1105 Input Parameters: 1106 + ts - The TS context obtained from TSCreate() 1107 - func - The function 1108 1109 Calling sequence of func: 1110 . func (TS ts); 1111 1112 Level: intermediate 1113 1114 .keywords: TS, timestep 1115 @*/ 1116 int TSSetPreStep(TS ts, int (*func)(TS)) 1117 { 1118 PetscFunctionBegin; 1119 PetscValidHeaderSpecific(ts, TS_COOKIE); 1120 ts->ops->prestep = func; 1121 PetscFunctionReturn(0); 1122 } 1123 1124 #undef __FUNC__ 1125 #define __FUNC__ "TSDefaultPreStep" 1126 /*@ 1127 TSDefaultPreStep - The default pre-stepping function which does nothing. 1128 1129 Collective on TS 1130 1131 Input Parameters: 1132 . ts - The TS context obtained from TSCreate() 1133 1134 Level: developer 1135 1136 .keywords: TS, timestep 1137 @*/ 1138 int TSDefaultPreStep(TS ts) 1139 { 1140 PetscFunctionBegin; 1141 PetscFunctionReturn(0); 1142 } 1143 1144 #undef __FUNC__ 1145 #define __FUNC__ "TSSetUpdate" 1146 /*@ 1147 TSSetUpdate - Sets the general-purpose update function called 1148 at the beginning of every time step. This function can change 1149 the time step. 1150 1151 Collective on TS 1152 1153 Input Parameters: 1154 + ts - The TS context obtained from TSCreate() 1155 - func - The function 1156 1157 Calling sequence of func: 1158 . func (TS ts, double t, double *dt); 1159 1160 + t - The current time 1161 - dt - The current time step 1162 1163 Level: intermediate 1164 1165 .keywords: TS, update, timestep 1166 @*/ 1167 int TSSetUpdate(TS ts, int (*func)(TS, double, double *)) 1168 { 1169 PetscFunctionBegin; 1170 PetscValidHeaderSpecific(ts, TS_COOKIE); 1171 ts->ops->update = func; 1172 PetscFunctionReturn(0); 1173 } 1174 1175 #undef __FUNC__ 1176 #define __FUNC__ "TSDefaultUpdate" 1177 /*@ 1178 TSDefaultUpdate - The default update function which does nothing. 1179 1180 Collective on TS 1181 1182 Input Parameters: 1183 + ts - The TS context obtained from TSCreate() 1184 - t - The current time 1185 1186 Output Parameters: 1187 . dt - The current time step 1188 1189 Level: developer 1190 1191 .keywords: TS, update, timestep 1192 @*/ 1193 int TSDefaultUpdate(TS ts, double t, double *dt) 1194 { 1195 PetscFunctionBegin; 1196 PetscFunctionReturn(0); 1197 } 1198 1199 #undef __FUNC__ 1200 #define __FUNC__ "TSSetPostStep" 1201 /*@ 1202 TSSetPostStep - Sets the general-purpose function 1203 called once at the end of time stepping. 1204 1205 Collective on TS 1206 1207 Input Parameters: 1208 + ts - The TS context obtained from TSCreate() 1209 - func - The function 1210 1211 Calling sequence of func: 1212 . func (TS ts); 1213 1214 Level: intermediate 1215 1216 .keywords: TS, timestep 1217 @*/ 1218 int TSSetPostStep(TS ts, int (*func)(TS)) 1219 { 1220 PetscFunctionBegin; 1221 PetscValidHeaderSpecific(ts, TS_COOKIE); 1222 ts->ops->poststep = func; 1223 PetscFunctionReturn(0); 1224 } 1225 1226 #undef __FUNC__ 1227 #define __FUNC__ "TSDefaultPostStep" 1228 /*@ 1229 TSDefaultPostStep - The default post-stepping function which does nothing. 1230 1231 Collective on TS 1232 1233 Input Parameters: 1234 . ts - The TS context obtained from TSCreate() 1235 1236 Level: developer 1237 1238 .keywords: TS, timestep 1239 @*/ 1240 int TSDefaultPostStep(TS ts) 1241 { 1242 PetscFunctionBegin; 1243 PetscFunctionReturn(0); 1244 } 1245 1246 /* ------------ Routines to set performance monitoring options ----------- */ 1247 1248 #undef __FUNCT__ 1249 #define __FUNCT__ "TSSetMonitor" 1250 /*@C 1251 TSSetMonitor - Sets an ADDITIONAL function that is to be used at every 1252 timestep to display the iteration's progress. 1253 1254 Collective on TS 1255 1256 Input Parameters: 1257 + ts - the TS context obtained from TSCreate() 1258 . func - monitoring routine 1259 . mctx - [optional] user-defined context for private data for the 1260 monitor routine (use PETSC_NULL if no context is desired) 1261 - monitordestroy - [optional] routine that frees monitor context 1262 (may be PETSC_NULL) 1263 1264 Calling sequence of func: 1265 $ int func(TS ts,int steps,PetscReal time,Vec x,void *mctx) 1266 1267 + ts - the TS context 1268 . steps - iteration number 1269 . time - current time 1270 . x - current iterate 1271 - mctx - [optional] monitoring context 1272 1273 Notes: 1274 This routine adds an additional monitor to the list of monitors that 1275 already has been loaded. 1276 1277 Level: intermediate 1278 1279 .keywords: TS, timestep, set, monitor 1280 1281 .seealso: TSDefaultMonitor(), TSClearMonitor() 1282 @*/ 1283 int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*)) 1284 { 1285 PetscFunctionBegin; 1286 PetscValidHeaderSpecific(ts,TS_COOKIE); 1287 if (ts->numbermonitors >= MAXTSMONITORS) { 1288 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1289 } 1290 ts->monitor[ts->numbermonitors] = monitor; 1291 ts->mdestroy[ts->numbermonitors] = mdestroy; 1292 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1293 PetscFunctionReturn(0); 1294 } 1295 1296 #undef __FUNCT__ 1297 #define __FUNCT__ "TSClearMonitor" 1298 /*@C 1299 TSClearMonitor - Clears all the monitors that have been set on a time-step object. 1300 1301 Collective on TS 1302 1303 Input Parameters: 1304 . ts - the TS context obtained from TSCreate() 1305 1306 Notes: 1307 There is no way to remove a single, specific monitor. 1308 1309 Level: intermediate 1310 1311 .keywords: TS, timestep, set, monitor 1312 1313 .seealso: TSDefaultMonitor(), TSSetMonitor() 1314 @*/ 1315 int TSClearMonitor(TS ts) 1316 { 1317 PetscFunctionBegin; 1318 PetscValidHeaderSpecific(ts,TS_COOKIE); 1319 ts->numbermonitors = 0; 1320 PetscFunctionReturn(0); 1321 } 1322 1323 #undef __FUNCT__ 1324 #define __FUNCT__ "TSDefaultMonitor" 1325 int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx) 1326 { 1327 int ierr; 1328 1329 PetscFunctionBegin; 1330 ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1331 PetscFunctionReturn(0); 1332 } 1333 1334 #undef __FUNCT__ 1335 #define __FUNCT__ "TSStep" 1336 /*@ 1337 TSStep - Steps the requested number of timesteps. 1338 1339 Collective on TS 1340 1341 Input Parameter: 1342 . ts - the TS context obtained from TSCreate() 1343 1344 Output Parameters: 1345 + steps - number of iterations until termination 1346 - ptime - time until termination 1347 1348 Level: beginner 1349 1350 .keywords: TS, timestep, solve 1351 1352 .seealso: TSCreate(), TSSetUp(), TSDestroy() 1353 @*/ 1354 int TSStep(TS ts,int *steps,PetscReal *ptime) 1355 { 1356 int ierr; 1357 PetscTruth flg; 1358 1359 PetscFunctionBegin; 1360 PetscValidHeaderSpecific(ts,TS_COOKIE); 1361 if (!ts->setupcalled) {ierr = TSSetUp(ts);CHKERRQ(ierr);} 1362 ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr); 1363 ierr = (*ts->ops->step)(ts,steps,ptime);CHKERRQ(ierr); 1364 ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr); 1365 ierr = PetscOptionsHasName(ts->prefix,"-ts_view",&flg);CHKERRQ(ierr); 1366 if (flg && !PetscPreLoadingOn) {ierr = TSView(ts,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 1367 PetscFunctionReturn(0); 1368 } 1369 1370 #undef __FUNCT__ 1371 #define __FUNCT__ "TSMonitor" 1372 /* 1373 Runs the user provided monitor routines, if they exists. 1374 */ 1375 int TSMonitor(TS ts,int step,PetscReal ptime,Vec x) 1376 { 1377 int i,ierr,n = ts->numbermonitors; 1378 1379 PetscFunctionBegin; 1380 for (i=0; i<n; i++) { 1381 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1382 } 1383 PetscFunctionReturn(0); 1384 } 1385 1386 /* ------------------------------------------------------------------------*/ 1387 1388 #undef __FUNCT__ 1389 #define __FUNCT__ "TSLGMonitorCreate" 1390 /*@C 1391 TSLGMonitorCreate - Creates a line graph context for use with 1392 TS to monitor convergence of preconditioned residual norms. 1393 1394 Collective on TS 1395 1396 Input Parameters: 1397 + host - the X display to open, or null for the local machine 1398 . label - the title to put in the title bar 1399 . x, y - the screen coordinates of the upper left coordinate of the window 1400 - m, n - the screen width and height in pixels 1401 1402 Output Parameter: 1403 . draw - the drawing context 1404 1405 Options Database Key: 1406 . -ts_xmonitor - automatically sets line graph monitor 1407 1408 Notes: 1409 Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1410 1411 Level: intermediate 1412 1413 .keywords: TS, monitor, line graph, residual, seealso 1414 1415 .seealso: TSLGMonitorDestroy(), TSSetMonitor() 1416 1417 @*/ 1418 int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw) 1419 { 1420 PetscDraw win; 1421 int ierr; 1422 1423 PetscFunctionBegin; 1424 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1425 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1426 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1427 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1428 1429 PetscLogObjectParent(*draw,win); 1430 PetscFunctionReturn(0); 1431 } 1432 1433 #undef __FUNCT__ 1434 #define __FUNCT__ "TSLGMonitor" 1435 int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx) 1436 { 1437 PetscDrawLG lg = (PetscDrawLG) monctx; 1438 PetscReal x,y = ptime; 1439 int ierr; 1440 1441 PetscFunctionBegin; 1442 if (!monctx) { 1443 MPI_Comm comm; 1444 PetscViewer viewer; 1445 1446 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1447 viewer = PETSC_VIEWER_DRAW_(comm); 1448 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1449 } 1450 1451 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1452 x = (PetscReal)n; 1453 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1454 if (n < 20 || (n % 5)) { 1455 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1456 } 1457 PetscFunctionReturn(0); 1458 } 1459 1460 #undef __FUNCT__ 1461 #define __FUNCT__ "TSLGMonitorDestroy" 1462 /*@C 1463 TSLGMonitorDestroy - Destroys a line graph context that was created 1464 with TSLGMonitorCreate(). 1465 1466 Collective on PetscDrawLG 1467 1468 Input Parameter: 1469 . draw - the drawing context 1470 1471 Level: intermediate 1472 1473 .keywords: TS, monitor, line graph, destroy 1474 1475 .seealso: TSLGMonitorCreate(), TSSetMonitor(), TSLGMonitor(); 1476 @*/ 1477 int TSLGMonitorDestroy(PetscDrawLG drawlg) 1478 { 1479 PetscDraw draw; 1480 int ierr; 1481 1482 PetscFunctionBegin; 1483 ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1484 ierr = PetscDrawDestroy(draw);CHKERRQ(ierr); 1485 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1486 PetscFunctionReturn(0); 1487 } 1488 1489 #undef __FUNCT__ 1490 #define __FUNCT__ "TSGetTime" 1491 /*@ 1492 TSGetTime - Gets the current time. 1493 1494 Not Collective 1495 1496 Input Parameter: 1497 . ts - the TS context obtained from TSCreate() 1498 1499 Output Parameter: 1500 . t - the current time 1501 1502 Contributed by: Matthew Knepley 1503 1504 Level: beginner 1505 1506 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1507 1508 .keywords: TS, get, time 1509 @*/ 1510 int TSGetTime(TS ts,PetscReal* t) 1511 { 1512 PetscFunctionBegin; 1513 PetscValidHeaderSpecific(ts,TS_COOKIE); 1514 *t = ts->ptime; 1515 PetscFunctionReturn(0); 1516 } 1517 1518 #undef __FUNCT__ 1519 #define __FUNCT__ "TSGetProblemType" 1520 /*@C 1521 TSGetProblemType - Returns the problem type of a TS (timestepper) context. 1522 1523 Not Collective 1524 1525 Input Parameter: 1526 . ts - The TS context obtained from TSCreate() 1527 1528 Output Parameter: 1529 . type - The problem type, TS_LINEAR or TS_NONLINEAR 1530 1531 Level: intermediate 1532 1533 Contributed by: Matthew Knepley 1534 1535 .keywords: ts, get, type 1536 1537 @*/ 1538 int TSGetProblemType(TS ts,TSProblemType *type) 1539 { 1540 PetscFunctionBegin; 1541 PetscValidHeaderSpecific(ts,TS_COOKIE); 1542 *type = ts->problem_type; 1543 PetscFunctionReturn(0); 1544 } 1545 1546 #undef __FUNCT__ 1547 #define __FUNCT__ "TSSetOptionsPrefix" 1548 /*@C 1549 TSSetOptionsPrefix - Sets the prefix used for searching for all 1550 TS options in the database. 1551 1552 Collective on TS 1553 1554 Input Parameter: 1555 + ts - The TS context 1556 - prefix - The prefix to prepend to all option names 1557 1558 Notes: 1559 A hyphen (-) must NOT be given at the beginning of the prefix name. 1560 The first character of all runtime options is AUTOMATICALLY the 1561 hyphen. 1562 1563 Contributed by: Matthew Knepley 1564 1565 Level: advanced 1566 1567 .keywords: TS, set, options, prefix, database 1568 1569 .seealso: TSSetFromOptions() 1570 1571 @*/ 1572 int TSSetOptionsPrefix(TS ts,char *prefix) 1573 { 1574 int ierr; 1575 1576 PetscFunctionBegin; 1577 PetscValidHeaderSpecific(ts,TS_COOKIE); 1578 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1579 switch(ts->problem_type) { 1580 case TS_NONLINEAR: 1581 ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1582 break; 1583 case TS_LINEAR: 1584 ierr = SLESSetOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr); 1585 break; 1586 } 1587 PetscFunctionReturn(0); 1588 } 1589 1590 1591 #undef __FUNCT__ 1592 #define __FUNCT__ "TSAppendOptionsPrefix" 1593 /*@C 1594 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1595 TS options in the database. 1596 1597 Collective on TS 1598 1599 Input Parameter: 1600 + ts - The TS context 1601 - prefix - The prefix to prepend to all option names 1602 1603 Notes: 1604 A hyphen (-) must NOT be given at the beginning of the prefix name. 1605 The first character of all runtime options is AUTOMATICALLY the 1606 hyphen. 1607 1608 Contributed by: Matthew Knepley 1609 1610 Level: advanced 1611 1612 .keywords: TS, append, options, prefix, database 1613 1614 .seealso: TSGetOptionsPrefix() 1615 1616 @*/ 1617 int TSAppendOptionsPrefix(TS ts,char *prefix) 1618 { 1619 int ierr; 1620 1621 PetscFunctionBegin; 1622 PetscValidHeaderSpecific(ts,TS_COOKIE); 1623 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1624 switch(ts->problem_type) { 1625 case TS_NONLINEAR: 1626 ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1627 break; 1628 case TS_LINEAR: 1629 ierr = SLESAppendOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr); 1630 break; 1631 } 1632 PetscFunctionReturn(0); 1633 } 1634 1635 #undef __FUNCT__ 1636 #define __FUNCT__ "TSGetOptionsPrefix" 1637 /*@C 1638 TSGetOptionsPrefix - Sets the prefix used for searching for all 1639 TS options in the database. 1640 1641 Not Collective 1642 1643 Input Parameter: 1644 . ts - The TS context 1645 1646 Output Parameter: 1647 . prefix - A pointer to the prefix string used 1648 1649 Contributed by: Matthew Knepley 1650 1651 Notes: On the fortran side, the user should pass in a string 'prifix' of 1652 sufficient length to hold the prefix. 1653 1654 Level: intermediate 1655 1656 .keywords: TS, get, options, prefix, database 1657 1658 .seealso: TSAppendOptionsPrefix() 1659 @*/ 1660 int TSGetOptionsPrefix(TS ts,char **prefix) 1661 { 1662 int ierr; 1663 1664 PetscFunctionBegin; 1665 PetscValidHeaderSpecific(ts,TS_COOKIE); 1666 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1667 PetscFunctionReturn(0); 1668 } 1669 1670 #undef __FUNCT__ 1671 #define __FUNCT__ "TSGetRHSMatrix" 1672 /*@C 1673 TSGetRHSMatrix - Returns the matrix A at the present timestep. 1674 1675 Not Collective, but parallel objects are returned if TS is parallel 1676 1677 Input Parameter: 1678 . ts - The TS context obtained from TSCreate() 1679 1680 Output Parameters: 1681 + A - The matrix A, where U_t = A(t) U 1682 . M - The preconditioner matrix, usually the same as A 1683 - ctx - User-defined context for matrix evaluation routine 1684 1685 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1686 1687 Contributed by: Matthew Knepley 1688 1689 Level: intermediate 1690 1691 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 1692 1693 .keywords: TS, timestep, get, matrix 1694 1695 @*/ 1696 int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx) 1697 { 1698 PetscFunctionBegin; 1699 PetscValidHeaderSpecific(ts,TS_COOKIE); 1700 if (A) *A = ts->A; 1701 if (M) *M = ts->B; 1702 if (ctx) *ctx = ts->jacP; 1703 PetscFunctionReturn(0); 1704 } 1705 1706 #undef __FUNCT__ 1707 #define __FUNCT__ "TSGetRHSJacobian" 1708 /*@C 1709 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1710 1711 Not Collective, but parallel objects are returned if TS is parallel 1712 1713 Input Parameter: 1714 . ts - The TS context obtained from TSCreate() 1715 1716 Output Parameters: 1717 + J - The Jacobian J of F, where U_t = F(U,t) 1718 . M - The preconditioner matrix, usually the same as J 1719 - ctx - User-defined context for Jacobian evaluation routine 1720 1721 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1722 1723 Contributed by: Matthew Knepley 1724 1725 Level: intermediate 1726 1727 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber() 1728 1729 .keywords: TS, timestep, get, matrix, Jacobian 1730 @*/ 1731 int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 1732 { 1733 int ierr; 1734 1735 PetscFunctionBegin; 1736 ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr); 1737 PetscFunctionReturn(0); 1738 } 1739 1740 /*MC 1741 TSRegisterDynamic - Adds a method to the timestepping solver package. 1742 1743 Synopsis: 1744 int TSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(TS)) 1745 1746 Not collective 1747 1748 Input Parameters: 1749 + name_solver - name of a new user-defined solver 1750 . path - path (either absolute or relative) the library containing this solver 1751 . name_create - name of routine to create method context 1752 - routine_create - routine to create method context 1753 1754 Notes: 1755 TSRegisterDynamic() may be called multiple times to add several user-defined solvers. 1756 1757 If dynamic libraries are used, then the fourth input argument (routine_create) 1758 is ignored. 1759 1760 Sample usage: 1761 .vb 1762 TSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 1763 "MySolverCreate",MySolverCreate); 1764 .ve 1765 1766 Then, your solver can be chosen with the procedural interface via 1767 $ TSSetType(ts,"my_solver") 1768 or at runtime via the option 1769 $ -ts_type my_solver 1770 1771 Level: advanced 1772 1773 Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT}, 1774 and others of the form ${any_environmental_variable} occuring in pathname will be 1775 replaced with appropriate values. 1776 1777 .keywords: TS, register 1778 1779 .seealso: TSRegisterAll(), TSRegisterDestroy() 1780 M*/ 1781 1782 #undef __FUNCT__ 1783 #define __FUNCT__ "TSRegister" 1784 int TSRegister(const char sname[], const char path[], const char name[], int (*function)(TS)) 1785 { 1786 char fullname[256]; 1787 int ierr; 1788 1789 PetscFunctionBegin; 1790 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 1791 ierr = PetscFListAdd(&TSList,sname,fullname,(void (*)())function);CHKERRQ(ierr); 1792 PetscFunctionReturn(0); 1793 } 1794 1795 #undef __FUNCT__ 1796 #define __FUNCT__ "TSVecViewMonitor" 1797 /*@C 1798 TSVecViewMonitor - Monitors progress of the TS solvers by calling 1799 VecView() for the solution at each timestep 1800 1801 Collective on TS 1802 1803 Input Parameters: 1804 + ts - the TS context 1805 . step - current time-step 1806 . ptime - current time 1807 - dummy - either a viewer or PETSC_NULL 1808 1809 Level: intermediate 1810 1811 .keywords: TS, vector, monitor, view 1812 1813 .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView() 1814 @*/ 1815 int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy) 1816 { 1817 int ierr; 1818 PetscViewer viewer = (PetscViewer) dummy; 1819 1820 PetscFunctionBegin; 1821 if (!viewer) { 1822 MPI_Comm comm; 1823 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1824 viewer = PETSC_VIEWER_DRAW_(comm); 1825 } 1826 ierr = VecView(x,viewer);CHKERRQ(ierr); 1827 PetscFunctionReturn(0); 1828 } 1829 1830 1831 1832