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 43 /*E 44 TSProblemType - Determines the type of problem this TS object is to be used to solve 45 46 Level: beginner 47 48 .seealso: TSCreate() 49 E*/ 50 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType; 51 52 /* Logging support */ 53 extern PetscClassId PETSCTS_DLLEXPORT TS_CLASSID; 54 55 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSInitializePackage(const char[]); 56 57 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate(MPI_Comm,TS*); 58 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS); 59 60 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS,TSProblemType); 61 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS,TSProblemType*); 62 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void*)); 63 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS); 64 65 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS,const char[]); 66 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS,const char[]); 67 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS,const char *[]); 68 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetFromOptions(TS); 69 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS); 70 71 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS,Vec); 72 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS,Vec*); 73 74 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS,PetscInt,PetscReal); 75 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS,PetscInt*,PetscReal*); 76 77 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 78 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*); 79 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS,PetscInt *,PetscReal*); 80 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS,Vec); 81 82 83 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS,PetscReal,PetscReal); 84 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS,PetscReal*); 85 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS,PetscReal*); 86 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS,PetscReal); 87 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS,PetscInt*); 88 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS,PetscReal); 89 90 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),void*); 91 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetMatrices(TS,Mat,PetscErrorCode (*)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat,PetscErrorCode (*)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure,void*); 92 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetMatrices(TS,Mat*,Mat*,void**); 93 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS,Mat,Mat,PetscErrorCode (*)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void*); 94 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS,Mat*,Mat*,void**); 95 96 typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 97 typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 98 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetIFunction(TS,TSIFunction,void*); 99 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 100 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 101 102 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 103 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 104 105 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS, PetscErrorCode (*)(TS)); 106 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS, PetscErrorCode (*)(TS)); 107 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS); 108 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS); 109 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPreStep(TS); 110 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPostStep(TS); 111 112 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 113 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultTimeStep(TS,PetscReal*,void*); 114 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPseudoComputeTimeStep(TS,PetscReal *); 115 116 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 117 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *); 118 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 119 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement(TS,PetscReal); 120 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt(TS); 121 122 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSPythonSetType(TS,const char[]); 123 124 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 125 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 126 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec); 127 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*); 128 129 /* Dynamic creation and loading functions */ 130 extern PetscFList TSList; 131 extern PetscBool TSRegisterAllCalled; 132 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetType(TS,const TSType*); 133 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetType(TS,const TSType); 134 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS)); 135 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSRegisterAll(const char[]); 136 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSRegisterDestroy(void); 137 138 /*MC 139 TSRegisterDynamic - Adds a creation method to the TS package. 140 141 Synopsis: 142 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 143 144 Not Collective 145 146 Input Parameters: 147 + name - The name of a new user-defined creation routine 148 . path - The path (either absolute or relative) of the library containing this routine 149 . func_name - The name of the creation routine 150 - create_func - The creation routine itself 151 152 Notes: 153 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 154 155 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 156 157 Sample usage: 158 .vb 159 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 160 .ve 161 162 Then, your ts type can be chosen with the procedural interface via 163 .vb 164 TS ts; 165 TSCreate(MPI_Comm, &ts); 166 TSSetType(ts, "my_ts") 167 .ve 168 or at runtime via the option 169 .vb 170 -ts_type my_ts 171 .ve 172 173 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 174 If your function is not being put into a shared library then use TSRegister() instead 175 176 Level: advanced 177 178 .keywords: TS, register 179 .seealso: TSRegisterAll(), TSRegisterDestroy() 180 M*/ 181 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 182 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 183 #else 184 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 185 #endif 186 187 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS,SNES*); 188 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS,KSP*); 189 190 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSView(TS,PetscViewer); 191 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSViewFromOptions(TS,const char[]); 192 193 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS,void *); 194 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS,void **); 195 196 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *); 197 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *); 198 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG); 199 200 /*S 201 TSGLAdapt - Abstract object that manages time-step adaptivity 202 203 Level: beginner 204 205 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 206 S*/ 207 typedef struct _p_TSGLAdapt *TSGLAdapt; 208 209 /*E 210 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 211 with an optional dynamic library name, for example 212 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 213 214 Level: beginner 215 216 .seealso: TSGLAdaptSetType(), TS 217 E*/ 218 #define TSGLAdaptType char* 219 #define TSGLADAPT_NONE "none" 220 #define TSGLADAPT_SIZE "size" 221 #define TSGLADAPT_BOTH "both" 222 223 /*MC 224 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 225 226 Synopsis: 227 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 228 229 Not Collective 230 231 Input Parameters: 232 + name_scheme - name of user-defined adaptivity scheme 233 . path - path (either absolute or relative) the library containing this scheme 234 . name_create - name of routine to create method context 235 - routine_create - routine to create method context 236 237 Notes: 238 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 239 240 If dynamic libraries are used, then the fourth input argument (routine_create) 241 is ignored. 242 243 Sample usage: 244 .vb 245 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 246 "MySchemeCreate",MySchemeCreate); 247 .ve 248 249 Then, your scheme can be chosen with the procedural interface via 250 $ TSGLAdaptSetType(ts,"my_scheme") 251 or at runtime via the option 252 $ -ts_adapt_type my_scheme 253 254 Level: advanced 255 256 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 257 and others of the form ${any_environmental_variable} occuring in pathname will be 258 replaced with appropriate values. 259 260 .keywords: TSGLAdapt, register 261 262 .seealso: TSGLAdaptRegisterAll() 263 M*/ 264 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 265 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 266 #else 267 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 268 #endif 269 270 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 271 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptRegisterAll(const char[]); 272 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptRegisterDestroy(void); 273 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptInitializePackage(const char[]); 274 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptFinalizePackage(void); 275 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 276 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType); 277 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 278 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 279 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptView(TSGLAdapt,PetscViewer); 280 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptSetFromOptions(TSGLAdapt); 281 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAdaptDestroy(TSGLAdapt); 282 283 /*E 284 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 285 with an optional dynamic library name, for example 286 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 287 288 Level: beginner 289 290 .seealso: TSGLSetAcceptType(), TS 291 E*/ 292 #define TSGLAcceptType char* 293 #define TSGLACCEPT_ALWAYS "always" 294 295 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 296 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 297 298 /*MC 299 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 300 301 Synopsis: 302 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 303 304 Not Collective 305 306 Input Parameters: 307 + name_scheme - name of user-defined acceptance scheme 308 . path - path (either absolute or relative) the library containing this scheme 309 . name_create - name of routine to create method context 310 - routine_create - routine to create method context 311 312 Notes: 313 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 314 315 If dynamic libraries are used, then the fourth input argument (routine_create) 316 is ignored. 317 318 Sample usage: 319 .vb 320 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 321 "MySchemeCreate",MySchemeCreate); 322 .ve 323 324 Then, your scheme can be chosen with the procedural interface via 325 $ TSGLSetAcceptType(ts,"my_scheme") 326 or at runtime via the option 327 $ -ts_gl_accept_type my_scheme 328 329 Level: advanced 330 331 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 332 and others of the form ${any_environmental_variable} occuring in pathname will be 333 replaced with appropriate values. 334 335 .keywords: TSGL, TSGLAcceptType, register 336 337 .seealso: TSGLRegisterAll() 338 M*/ 339 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 340 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 341 #else 342 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 343 #endif 344 345 /*E 346 TSGLType - family of time integration method within the General Linear class 347 348 Level: beginner 349 350 .seealso: TSGLSetType(), TSGLRegister() 351 E*/ 352 #define TSGLType char* 353 #define TSGL_IRKS "irks" 354 355 /*MC 356 TSGLRegisterDynamic - adds a TSGL implementation 357 358 Synopsis: 359 PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 360 361 Not Collective 362 363 Input Parameters: 364 + name_scheme - name of user-defined general linear scheme 365 . path - path (either absolute or relative) the library containing this scheme 366 . name_create - name of routine to create method context 367 - routine_create - routine to create method context 368 369 Notes: 370 TSGLRegisterDynamic() may be called multiple times to add several user-defined families. 371 372 If dynamic libraries are used, then the fourth input argument (routine_create) 373 is ignored. 374 375 Sample usage: 376 .vb 377 TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 378 "MySchemeCreate",MySchemeCreate); 379 .ve 380 381 Then, your scheme can be chosen with the procedural interface via 382 $ TSGLSetType(ts,"my_scheme") 383 or at runtime via the option 384 $ -ts_gl_type my_scheme 385 386 Level: advanced 387 388 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 389 and others of the form ${any_environmental_variable} occuring in pathname will be 390 replaced with appropriate values. 391 392 .keywords: TSGL, register 393 394 .seealso: TSGLRegisterAll() 395 M*/ 396 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 397 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 398 #else 399 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 400 #endif 401 402 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 403 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLRegisterAll(const char[]); 404 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLRegisterDestroy(void); 405 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLInitializePackage(const char[]); 406 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLFinalizePackage(void); 407 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLSetType(TS,const TSGLType); 408 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLGetAdapt(TS,TSGLAdapt*); 409 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGLSetAcceptType(TS,const TSGLAcceptType); 410 411 /* 412 PETSc interface to Sundials 413 */ 414 #ifdef PETSC_HAVE_SUNDIALS 415 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 416 extern const char *TSSundialsLmmTypes[]; 417 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 418 extern const char *TSSundialsGramSchmidtTypes[]; 419 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS,TSSundialsLmmType); 420 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS,PC*); 421 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS,PetscReal,PetscReal); 422 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMinTimeStep(TS,PetscReal); 423 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMaxTimeStep(TS,PetscReal); 424 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 425 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 426 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS,PetscInt); 427 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS,PetscReal); 428 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS,PetscBool ); 429 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps(TS,PetscBool ); 430 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 431 #endif 432 433 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSRKSetTolerance(TS,PetscReal); 434 435 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSThetaSetTheta(TS,PetscReal); 436 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSThetaGetTheta(TS,PetscReal*); 437 438 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSAlphaSetAccept(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 439 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSAlphaSetRadius(TS,PetscReal); 440 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 441 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 442 443 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSSetDM(TS,DM); 444 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSGetDM(TS,DM*); 445 446 EXTERN PetscErrorCode PETSCTS_DLLEXPORT SNESTSFormFunction(SNES,Vec,Vec,void*); 447 EXTERN PetscErrorCode PETSCTS_DLLEXPORT SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 448 449 PETSC_EXTERN_CXX_END 450 #endif 451