1 /* 2 User interface for the timestepping package. This package 3 is for use in solving time-dependent PDEs. 4 */ 5 #if !defined(__PETSCTS_H) 6 #define __PETSCTS_H 7 #include "petscsnes.h" 8 PETSC_EXTERN_CXX_BEGIN 9 10 /*S 11 TS - Abstract PETSc object that manages all time-steppers (ODE integrators) 12 13 Level: beginner 14 15 Concepts: ODE solvers 16 17 .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC 18 S*/ 19 typedef struct _p_TS* TS; 20 21 /*E 22 TSType - String with the name of a PETSc TS method or the creation function 23 with an optional dynamic library name, for example 24 http://www.mcs.anl.gov/petsc/lib.a:mytscreate() 25 26 Level: beginner 27 28 .seealso: TSSetType(), TS 29 E*/ 30 #define TSType char* 31 #define TSEULER "euler" 32 #define TSBEULER "beuler" 33 #define TSPSEUDO "pseudo" 34 #define TSCN "cn" 35 #define TSSUNDIALS "sundials" 36 #define TSRK "rk" 37 #define TSPYTHON "python" 38 #define TSTHETA "theta" 39 #define TSALPHA "alpha" 40 #define TSGL "gl" 41 #define TSSSP "ssp" 42 #define TSARKIMEX "arkimex" 43 #if defined(PETSC_HAVE_ROSW) 44 #define TSROSW "rosw" 45 #endif 46 47 /*E 48 TSProblemType - Determines the type of problem this TS object is to be used to solve 49 50 Level: beginner 51 52 .seealso: TSCreate() 53 E*/ 54 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType; 55 56 typedef enum { 57 TS_CONVERGED_ITERATING = 0, 58 TS_CONVERGED_TIME = 1, 59 TS_CONVERGED_ITS = 2, 60 TS_DIVERGED_NONLINEAR_SOLVE = -1, 61 TS_DIVERGED_STEP_REJECTED = -2 62 } TSConvergedReason; 63 extern const char *const*TSConvergedReasons; 64 65 /* Logging support */ 66 extern PetscClassId TS_CLASSID; 67 68 extern PetscErrorCode TSInitializePackage(const char[]); 69 70 extern PetscErrorCode TSCreate(MPI_Comm,TS*); 71 extern PetscErrorCode TSDestroy(TS*); 72 73 extern PetscErrorCode TSSetProblemType(TS,TSProblemType); 74 extern PetscErrorCode TSGetProblemType(TS,TSProblemType*); 75 extern PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 76 extern PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 77 extern PetscErrorCode TSMonitorCancel(TS); 78 79 extern PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 80 extern PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 81 extern PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 82 extern PetscErrorCode TSSetFromOptions(TS); 83 extern PetscErrorCode TSSetUp(TS); 84 extern PetscErrorCode TSReset(TS); 85 86 extern PetscErrorCode TSSetSolution(TS,Vec); 87 extern PetscErrorCode TSGetSolution(TS,Vec*); 88 89 extern PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 90 extern PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 91 extern PetscErrorCode TSSetExactFinalTime(TS,PetscBool); 92 93 extern PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 94 extern PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*); 95 extern PetscErrorCode TSMonitorSolutionCreate(TS,PetscViewer,void**); 96 extern PetscErrorCode TSMonitorSolutionDestroy(void**); 97 extern PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*); 98 99 extern PetscErrorCode TSStep(TS); 100 extern PetscErrorCode TSSolve(TS,Vec,PetscReal*); 101 extern PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 102 103 extern PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 104 extern PetscErrorCode TSGetTimeStep(TS,PetscReal*); 105 extern PetscErrorCode TSGetTime(TS,PetscReal*); 106 extern PetscErrorCode TSSetTime(TS,PetscReal); 107 extern PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 108 extern PetscErrorCode TSSetTimeStep(TS,PetscReal); 109 110 typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 111 typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 112 extern PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 113 extern PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 114 extern PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 115 extern PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 116 117 typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 118 typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 119 extern PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 120 extern PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 121 extern PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 122 extern PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 123 124 extern PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 125 extern PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 126 extern PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 127 extern PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 128 extern PetscErrorCode TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 129 extern PetscErrorCode TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 130 131 extern PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 132 extern PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 133 extern PetscErrorCode TSPreStep(TS); 134 extern PetscErrorCode TSPostStep(TS); 135 extern PetscErrorCode TSSetRetainStages(TS,PetscBool); 136 extern PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 137 138 extern PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 139 extern PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*); 140 extern PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 141 142 extern PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 143 extern PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *); 144 extern PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 145 extern PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 146 extern PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 147 148 extern PetscErrorCode TSPythonSetType(TS,const char[]); 149 150 extern PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 151 extern PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 152 extern PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 153 extern PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool); 154 155 extern PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 156 157 /* Dynamic creation and loading functions */ 158 extern PetscFList TSList; 159 extern PetscBool TSRegisterAllCalled; 160 extern PetscErrorCode TSGetType(TS,const TSType*); 161 extern PetscErrorCode TSSetType(TS,const TSType); 162 extern PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS)); 163 extern PetscErrorCode TSRegisterAll(const char[]); 164 extern PetscErrorCode TSRegisterDestroy(void); 165 166 /*MC 167 TSRegisterDynamic - Adds a creation method to the TS package. 168 169 Synopsis: 170 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 171 172 Not Collective 173 174 Input Parameters: 175 + name - The name of a new user-defined creation routine 176 . path - The path (either absolute or relative) of the library containing this routine 177 . func_name - The name of the creation routine 178 - create_func - The creation routine itself 179 180 Notes: 181 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 182 183 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 184 185 Sample usage: 186 .vb 187 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 188 .ve 189 190 Then, your ts type can be chosen with the procedural interface via 191 .vb 192 TS ts; 193 TSCreate(MPI_Comm, &ts); 194 TSSetType(ts, "my_ts") 195 .ve 196 or at runtime via the option 197 .vb 198 -ts_type my_ts 199 .ve 200 201 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 202 If your function is not being put into a shared library then use TSRegister() instead 203 204 Level: advanced 205 206 .keywords: TS, register 207 .seealso: TSRegisterAll(), TSRegisterDestroy() 208 M*/ 209 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 210 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 211 #else 212 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 213 #endif 214 215 extern PetscErrorCode TSGetSNES(TS,SNES*); 216 extern PetscErrorCode TSGetKSP(TS,KSP*); 217 218 extern PetscErrorCode TSView(TS,PetscViewer); 219 220 extern PetscErrorCode TSSetApplicationContext(TS,void *); 221 extern PetscErrorCode TSGetApplicationContext(TS,void *); 222 223 extern PetscErrorCode TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *); 224 extern PetscErrorCode TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *); 225 extern PetscErrorCode TSMonitorLGDestroy(PetscDrawLG*); 226 227 /*E 228 TSSSPType - string with the name of TSSSP scheme. 229 230 Level: beginner 231 232 .seealso: TSSSPSetType(), TS 233 E*/ 234 #define TSSSPType char* 235 #define TSSSPRKS2 "rks2" 236 #define TSSSPRKS3 "rks3" 237 #define TSSSPRK104 "rk104" 238 239 extern PetscErrorCode TSSSPSetType(TS,const TSSSPType); 240 extern PetscErrorCode TSSSPGetType(TS,const TSSSPType*); 241 extern PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 242 extern PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 243 244 /*S 245 TSGLAdapt - Abstract object that manages time-step adaptivity 246 247 Level: beginner 248 249 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 250 S*/ 251 typedef struct _p_TSGLAdapt *TSGLAdapt; 252 253 /*E 254 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 255 with an optional dynamic library name, for example 256 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 257 258 Level: beginner 259 260 .seealso: TSGLAdaptSetType(), TS 261 E*/ 262 #define TSGLAdaptType char* 263 #define TSGLADAPT_NONE "none" 264 #define TSGLADAPT_SIZE "size" 265 #define TSGLADAPT_BOTH "both" 266 267 /*MC 268 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 269 270 Synopsis: 271 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 272 273 Not Collective 274 275 Input Parameters: 276 + name_scheme - name of user-defined adaptivity scheme 277 . path - path (either absolute or relative) the library containing this scheme 278 . name_create - name of routine to create method context 279 - routine_create - routine to create method context 280 281 Notes: 282 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 283 284 If dynamic libraries are used, then the fourth input argument (routine_create) 285 is ignored. 286 287 Sample usage: 288 .vb 289 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 290 "MySchemeCreate",MySchemeCreate); 291 .ve 292 293 Then, your scheme can be chosen with the procedural interface via 294 $ TSGLAdaptSetType(ts,"my_scheme") 295 or at runtime via the option 296 $ -ts_adapt_type my_scheme 297 298 Level: advanced 299 300 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 301 and others of the form ${any_environmental_variable} occuring in pathname will be 302 replaced with appropriate values. 303 304 .keywords: TSGLAdapt, register 305 306 .seealso: TSGLAdaptRegisterAll() 307 M*/ 308 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 309 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 310 #else 311 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 312 #endif 313 314 extern PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 315 extern PetscErrorCode TSGLAdaptRegisterAll(const char[]); 316 extern PetscErrorCode TSGLAdaptRegisterDestroy(void); 317 extern PetscErrorCode TSGLAdaptInitializePackage(const char[]); 318 extern PetscErrorCode TSGLAdaptFinalizePackage(void); 319 extern PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 320 extern PetscErrorCode TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType); 321 extern PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 322 extern PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 323 extern PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 324 extern PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 325 extern PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 326 327 /*E 328 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 329 with an optional dynamic library name, for example 330 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 331 332 Level: beginner 333 334 .seealso: TSGLSetAcceptType(), TS 335 E*/ 336 #define TSGLAcceptType char* 337 #define TSGLACCEPT_ALWAYS "always" 338 339 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 340 extern PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 341 342 /*MC 343 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 344 345 Synopsis: 346 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 347 348 Not Collective 349 350 Input Parameters: 351 + name_scheme - name of user-defined acceptance scheme 352 . path - path (either absolute or relative) the library containing this scheme 353 . name_create - name of routine to create method context 354 - routine_create - routine to create method context 355 356 Notes: 357 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 358 359 If dynamic libraries are used, then the fourth input argument (routine_create) 360 is ignored. 361 362 Sample usage: 363 .vb 364 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 365 "MySchemeCreate",MySchemeCreate); 366 .ve 367 368 Then, your scheme can be chosen with the procedural interface via 369 $ TSGLSetAcceptType(ts,"my_scheme") 370 or at runtime via the option 371 $ -ts_gl_accept_type my_scheme 372 373 Level: advanced 374 375 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 376 and others of the form ${any_environmental_variable} occuring in pathname will be 377 replaced with appropriate values. 378 379 .keywords: TSGL, TSGLAcceptType, register 380 381 .seealso: TSGLRegisterAll() 382 M*/ 383 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 384 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 385 #else 386 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 387 #endif 388 389 /*E 390 TSGLType - family of time integration method within the General Linear class 391 392 Level: beginner 393 394 .seealso: TSGLSetType(), TSGLRegister() 395 E*/ 396 #define TSGLType char* 397 #define TSGL_IRKS "irks" 398 399 /*MC 400 TSGLRegisterDynamic - adds a TSGL implementation 401 402 Synopsis: 403 PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 404 405 Not Collective 406 407 Input Parameters: 408 + name_scheme - name of user-defined general linear scheme 409 . path - path (either absolute or relative) the library containing this scheme 410 . name_create - name of routine to create method context 411 - routine_create - routine to create method context 412 413 Notes: 414 TSGLRegisterDynamic() may be called multiple times to add several user-defined families. 415 416 If dynamic libraries are used, then the fourth input argument (routine_create) 417 is ignored. 418 419 Sample usage: 420 .vb 421 TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 422 "MySchemeCreate",MySchemeCreate); 423 .ve 424 425 Then, your scheme can be chosen with the procedural interface via 426 $ TSGLSetType(ts,"my_scheme") 427 or at runtime via the option 428 $ -ts_gl_type my_scheme 429 430 Level: advanced 431 432 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 433 and others of the form ${any_environmental_variable} occuring in pathname will be 434 replaced with appropriate values. 435 436 .keywords: TSGL, register 437 438 .seealso: TSGLRegisterAll() 439 M*/ 440 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 441 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 442 #else 443 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 444 #endif 445 446 extern PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 447 extern PetscErrorCode TSGLRegisterAll(const char[]); 448 extern PetscErrorCode TSGLRegisterDestroy(void); 449 extern PetscErrorCode TSGLInitializePackage(const char[]); 450 extern PetscErrorCode TSGLFinalizePackage(void); 451 extern PetscErrorCode TSGLSetType(TS,const TSGLType); 452 extern PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 453 extern PetscErrorCode TSGLSetAcceptType(TS,const TSGLAcceptType); 454 455 #define TSARKIMEXType char* 456 #define TSARKIMEX2D "2d" 457 #define TSARKIMEX2E "2e" 458 #define TSARKIMEX3 "3" 459 #define TSARKIMEX4 "4" 460 #define TSARKIMEX5 "5" 461 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*); 462 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType); 463 extern PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 464 extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]); 465 extern PetscErrorCode TSARKIMEXFinalizePackage(void); 466 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]); 467 extern PetscErrorCode TSARKIMEXRegisterDestroy(void); 468 extern PetscErrorCode TSARKIMEXRegisterAll(void); 469 470 #if defined(PETSC_HAVE_ROSW) 471 #define TSRosWType char* 472 #define TSROSW2D "2d" 473 #define TSROSW2E "2e" 474 #define TSROSW3 "3" 475 #define TSROSW4 "4" 476 #define TSROSW5 "5" 477 extern PetscErrorCode TSRosWGetType(TS ts,const TSRosWType*); 478 extern PetscErrorCode TSRosWSetType(TS ts,const TSRosWType); 479 extern PetscErrorCode TSRosWSetFullyImplicit(TS,PetscBool); 480 extern PetscErrorCode TSRosWRegister(const TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]); 481 extern PetscErrorCode TSRosWFinalizePackage(void); 482 extern PetscErrorCode TSRosWInitializePackage(const char path[]); 483 extern PetscErrorCode TSRosWRegisterDestroy(void); 484 extern PetscErrorCode TSRosWRegisterAll(void); 485 #endif 486 487 /* 488 PETSc interface to Sundials 489 */ 490 #ifdef PETSC_HAVE_SUNDIALS 491 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 492 extern const char *TSSundialsLmmTypes[]; 493 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 494 extern const char *TSSundialsGramSchmidtTypes[]; 495 extern PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 496 extern PetscErrorCode TSSundialsGetPC(TS,PC*); 497 extern PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 498 extern PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 499 extern PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 500 extern PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 501 extern PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 502 extern PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 503 extern PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 504 extern PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 505 extern PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 506 extern PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 507 #endif 508 509 extern PetscErrorCode TSRKSetTolerance(TS,PetscReal); 510 511 extern PetscErrorCode TSThetaSetTheta(TS,PetscReal); 512 extern PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 513 extern PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 514 extern PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 515 516 extern PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 517 extern PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 518 extern PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 519 extern PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 520 extern PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 521 522 extern PetscErrorCode TSSetDM(TS,DM); 523 extern PetscErrorCode TSGetDM(TS,DM*); 524 525 extern PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 526 extern PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 527 528 PETSC_EXTERN_CXX_END 529 #endif 530