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