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