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