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