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