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