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