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