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