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