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