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