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