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