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 9 /*S 10 TS - Abstract PETSc object that manages all time-steppers (ODE integrators) 11 12 Level: beginner 13 14 Concepts: ODE solvers 15 16 .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC 17 S*/ 18 typedef struct _p_TS* TS; 19 20 /*J 21 TSType - String with the name of a PETSc TS method or the creation function 22 with an optional dynamic library name, for example 23 http://www.mcs.anl.gov/petsc/lib.a:mytscreate() 24 25 Level: beginner 26 27 .seealso: TSSetType(), TS 28 J*/ 29 typedef const char* TSType; 30 #define TSEULER "euler" 31 #define TSBEULER "beuler" 32 #define TSPSEUDO "pseudo" 33 #define TSCN "cn" 34 #define TSSUNDIALS "sundials" 35 #define TSRK "rk" 36 #define TSPYTHON "python" 37 #define TSTHETA "theta" 38 #define TSALPHA "alpha" 39 #define TSGL "gl" 40 #define TSSSP "ssp" 41 #define TSARKIMEX "arkimex" 42 #define TSROSW "rosw" 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 /*E 54 TSEquationType - type of TS problem that is solved 55 56 Level: beginner 57 58 Developer Notes: this must match finclude/petscts.h 59 60 Supported types are: 61 TS_EQ_UNSPECIFIED (default) 62 TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0 63 TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0 64 65 .seealso: TSGetEquationType(), TSSetEquationType() 66 E*/ 67 typedef enum { 68 TS_EQ_UNSPECIFIED = -1, 69 TS_EQ_EXPLICIT = 0, 70 TS_EQ_ODE_EXPLICIT = 1, 71 TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100, 72 TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200, 73 TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300, 74 TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500, 75 TS_EQ_IMPLICIT = 1000, 76 TS_EQ_ODE_IMPLICIT = 1001, 77 TS_EQ_DAE_IMPLICIT_INDEX1 = 1100, 78 TS_EQ_DAE_IMPLICIT_INDEX2 = 1200, 79 TS_EQ_DAE_IMPLICIT_INDEX3 = 1300, 80 TS_EQ_DAE_IMPLICIT_INDEXHI = 1500, 81 } TSEquationType; 82 PETSC_EXTERN const char *const*TSEquationTypes; 83 84 /*E 85 TSConvergedReason - reason a TS method has converged or not 86 87 Level: beginner 88 89 Developer Notes: this must match finclude/petscts.h 90 91 Each reason has its own manual page. 92 93 .seealso: TSGetConvergedReason() 94 E*/ 95 typedef enum { 96 TS_CONVERGED_ITERATING = 0, 97 TS_CONVERGED_TIME = 1, 98 TS_CONVERGED_ITS = 2, 99 TS_CONVERGED_USER = 3, 100 TS_DIVERGED_NONLINEAR_SOLVE = -1, 101 TS_DIVERGED_STEP_REJECTED = -2 102 } TSConvergedReason; 103 PETSC_EXTERN const char *const*TSConvergedReasons; 104 /*MC 105 TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve() 106 107 Level: beginner 108 109 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt() 110 M*/ 111 112 /*MC 113 TS_CONVERGED_TIME - the final time was reached 114 115 Level: beginner 116 117 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt(), TSSetDuration(), TSGetSolveTime() 118 M*/ 119 120 /*MC 121 TS_CONVERGED_ITS - the maximum number of iterations was reached prior to the final time 122 123 Level: beginner 124 125 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt(), TSSetDuration() 126 M*/ 127 /*MC 128 TS_CONVERGED_USER - user requested termination 129 130 Level: beginner 131 132 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 133 M*/ 134 135 /*MC 136 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 137 138 Level: beginner 139 140 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt(), TSGetSNES(), SNESGetConvergedReason() 141 M*/ 142 143 /*MC 144 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 145 146 Level: beginner 147 148 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt() 149 M*/ 150 151 /*E 152 TSExactFinalTimeOption - option for handling of final time step 153 154 Level: beginner 155 156 Developer Notes: this must match finclude/petscts.h 157 158 $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded 159 $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time 160 $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time 161 .seealso: TSGetConvergedReason() 162 E*/ 163 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption; 164 PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 165 166 167 /* Logging support */ 168 PETSC_EXTERN PetscClassId TS_CLASSID; 169 PETSC_EXTERN PetscClassId DMTS_CLASSID; 170 171 PETSC_EXTERN PetscErrorCode TSInitializePackage(const char[]); 172 173 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*); 174 PETSC_EXTERN PetscErrorCode TSDestroy(TS*); 175 176 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType); 177 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*); 178 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 179 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 180 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 181 182 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 183 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 184 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 185 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 186 PETSC_EXTERN PetscErrorCode TSSetUp(TS); 187 PETSC_EXTERN PetscErrorCode TSReset(TS); 188 189 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec); 190 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*); 191 192 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 193 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 194 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption); 195 196 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 197 198 typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx; 199 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *); 200 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*); 201 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*); 202 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*); 203 204 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*); 205 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 206 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 207 208 PETSC_EXTERN PetscErrorCode TSStep(TS); 209 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 210 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec); 211 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*); 212 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType); 213 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 214 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason); 215 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*); 216 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*); 217 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*); 218 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*); 219 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt); 220 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*); 221 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt); 222 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool); 223 224 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 225 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*); 226 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*); 227 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal); 228 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 229 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal); 230 231 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 232 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 233 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 234 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 235 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 236 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 237 238 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*); 239 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*); 240 241 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 242 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 243 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 244 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 245 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 246 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 247 248 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 249 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 250 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 251 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 252 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec); 253 254 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 255 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal)); 256 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 257 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 258 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal); 259 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 260 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool); 261 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 262 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 263 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*); 264 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*); 265 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 266 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*); 267 268 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 269 PETSC_EXTERN PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*); 270 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 271 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 272 273 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 274 PETSC_EXTERN PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *); 275 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 276 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 277 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 278 279 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]); 280 281 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 282 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 283 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 284 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool); 285 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*); 286 287 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 288 289 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*); 290 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**); 291 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*); 292 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**); 293 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*); 294 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**); 295 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*); 296 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**); 297 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*); 298 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**); 299 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 300 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 301 302 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 303 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*); 304 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 305 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*); 306 307 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *); 308 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*),void *); 309 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *); 310 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*),void *); 311 312 typedef struct { 313 Vec ray; 314 VecScatter scatter; 315 PetscViewer viewer; 316 } TSMonitorDMDARayCtx; 317 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**); 318 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*); 319 320 321 /* Dynamic creation and loading functions */ 322 PETSC_EXTERN PetscFList TSList; 323 PETSC_EXTERN PetscBool TSRegisterAllCalled; 324 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*); 325 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType); 326 PETSC_EXTERN PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS)); 327 PETSC_EXTERN PetscErrorCode TSRegisterAll(const char[]); 328 PETSC_EXTERN PetscErrorCode TSRegisterDestroy(void); 329 330 /*MC 331 TSRegisterDynamic - Adds a creation method to the TS package. 332 333 Synopsis: 334 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 335 336 Not Collective 337 338 Input Parameters: 339 + name - The name of a new user-defined creation routine 340 . path - The path (either absolute or relative) of the library containing this routine 341 . func_name - The name of the creation routine 342 - create_func - The creation routine itself 343 344 Notes: 345 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 346 347 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 348 349 Sample usage: 350 .vb 351 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 352 .ve 353 354 Then, your ts type can be chosen with the procedural interface via 355 .vb 356 TS ts; 357 TSCreate(MPI_Comm, &ts); 358 TSSetType(ts, "my_ts") 359 .ve 360 or at runtime via the option 361 .vb 362 -ts_type my_ts 363 .ve 364 365 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 366 If your function is not being put into a shared library then use TSRegister() instead 367 368 Level: advanced 369 370 .keywords: TS, register 371 .seealso: TSRegisterAll(), TSRegisterDestroy() 372 M*/ 373 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 374 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 375 #else 376 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 377 #endif 378 379 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 380 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES); 381 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 382 383 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 384 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer); 385 386 #define TS_FILE_CLASSID 1211225 387 388 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 389 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 390 391 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 392 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 393 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 394 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 395 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 396 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 397 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 398 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 399 400 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 401 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 402 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 403 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 404 405 /*J 406 TSSSPType - string with the name of TSSSP scheme. 407 408 Level: beginner 409 410 .seealso: TSSSPSetType(), TS 411 J*/ 412 typedef const char* TSSSPType; 413 #define TSSSPRKS2 "rks2" 414 #define TSSSPRKS3 "rks3" 415 #define TSSSPRK104 "rk104" 416 417 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 418 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 419 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 420 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 421 422 /*S 423 TSAdapt - Abstract object that manages time-step adaptivity 424 425 Level: beginner 426 427 .seealso: TS, TSAdaptCreate(), TSAdaptType 428 S*/ 429 typedef struct _p_TSAdapt *TSAdapt; 430 431 /*E 432 TSAdaptType - String with the name of TSAdapt scheme or the creation function 433 with an optional dynamic library name, for example 434 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 435 436 Level: beginner 437 438 .seealso: TSAdaptSetType(), TS 439 E*/ 440 typedef const char *TSAdaptType; 441 #define TSADAPTBASIC "basic" 442 #define TSADAPTNONE "none" 443 #define TSADAPTCFL "cfl" 444 445 /*MC 446 TSAdaptRegisterDynamic - adds a TSAdapt implementation 447 448 Synopsis: 449 PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 450 451 Not Collective 452 453 Input Parameters: 454 + name_scheme - name of user-defined adaptivity scheme 455 . path - path (either absolute or relative) the library containing this scheme 456 . name_create - name of routine to create method context 457 - routine_create - routine to create method context 458 459 Notes: 460 TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 461 462 If dynamic libraries are used, then the fourth input argument (routine_create) 463 is ignored. 464 465 Sample usage: 466 .vb 467 TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 468 "MySchemeCreate",MySchemeCreate); 469 .ve 470 471 Then, your scheme can be chosen with the procedural interface via 472 $ TSAdaptSetType(ts,"my_scheme") 473 or at runtime via the option 474 $ -ts_adapt_type my_scheme 475 476 Level: advanced 477 478 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 479 and others of the form ${any_environmental_variable} occuring in pathname will be 480 replaced with appropriate values. 481 482 .keywords: TSAdapt, register 483 484 .seealso: TSAdaptRegisterAll() 485 M*/ 486 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 487 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,0) 488 #else 489 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,d) 490 #endif 491 492 PETSC_EXTERN PetscErrorCode TSGetTSAdapt(TS,TSAdapt*); 493 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt)); 494 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(const char[]); 495 PETSC_EXTERN PetscErrorCode TSAdaptRegisterDestroy(void); 496 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(const char[]); 497 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 498 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 499 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 500 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 501 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 502 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 503 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 504 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 505 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*); 506 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 507 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer); 508 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt); 509 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 510 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 511 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 512 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*)); 513 514 /*S 515 TSGLAdapt - Abstract object that manages time-step adaptivity 516 517 Level: beginner 518 519 Developer Notes: 520 This functionality should be replaced by the TSAdapt. 521 522 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 523 S*/ 524 typedef struct _p_TSGLAdapt *TSGLAdapt; 525 526 /*J 527 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 528 with an optional dynamic library name, for example 529 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 530 531 Level: beginner 532 533 .seealso: TSGLAdaptSetType(), TS 534 J*/ 535 typedef const char *TSGLAdaptType; 536 #define TSGLADAPT_NONE "none" 537 #define TSGLADAPT_SIZE "size" 538 #define TSGLADAPT_BOTH "both" 539 540 /*MC 541 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 542 543 Synopsis: 544 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 545 546 Not Collective 547 548 Input Parameters: 549 + name_scheme - name of user-defined adaptivity scheme 550 . path - path (either absolute or relative) the library containing this scheme 551 . name_create - name of routine to create method context 552 - routine_create - routine to create method context 553 554 Notes: 555 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 556 557 If dynamic libraries are used, then the fourth input argument (routine_create) 558 is ignored. 559 560 Sample usage: 561 .vb 562 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 563 "MySchemeCreate",MySchemeCreate); 564 .ve 565 566 Then, your scheme can be chosen with the procedural interface via 567 $ TSGLAdaptSetType(ts,"my_scheme") 568 or at runtime via the option 569 $ -ts_adapt_type my_scheme 570 571 Level: advanced 572 573 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 574 and others of the form ${any_environmental_variable} occuring in pathname will be 575 replaced with appropriate values. 576 577 .keywords: TSGLAdapt, register 578 579 .seealso: TSGLAdaptRegisterAll() 580 M*/ 581 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 582 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 583 #else 584 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 585 #endif 586 587 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 588 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(const char[]); 589 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterDestroy(void); 590 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(const char[]); 591 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void); 592 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 593 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType); 594 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 595 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 596 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 597 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 598 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 599 600 /*J 601 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 602 with an optional dynamic library name, for example 603 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 604 605 Level: beginner 606 607 .seealso: TSGLSetAcceptType(), TS 608 J*/ 609 typedef const char *TSGLAcceptType; 610 #define TSGLACCEPT_ALWAYS "always" 611 612 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 613 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 614 615 /*MC 616 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 617 618 Synopsis: 619 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 620 621 Not Collective 622 623 Input Parameters: 624 + name_scheme - name of user-defined acceptance scheme 625 . path - path (either absolute or relative) the library containing this scheme 626 . name_create - name of routine to create method context 627 - routine_create - routine to create method context 628 629 Notes: 630 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 631 632 If dynamic libraries are used, then the fourth input argument (routine_create) 633 is ignored. 634 635 Sample usage: 636 .vb 637 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 638 "MySchemeCreate",MySchemeCreate); 639 .ve 640 641 Then, your scheme can be chosen with the procedural interface via 642 $ TSGLSetAcceptType(ts,"my_scheme") 643 or at runtime via the option 644 $ -ts_gl_accept_type my_scheme 645 646 Level: advanced 647 648 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 649 and others of the form ${any_environmental_variable} occuring in pathname will be 650 replaced with appropriate values. 651 652 .keywords: TSGL, TSGLAcceptType, register 653 654 .seealso: TSGLRegisterAll() 655 M*/ 656 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 657 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 658 #else 659 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 660 #endif 661 662 /*J 663 TSGLType - family of time integration method within the General Linear class 664 665 Level: beginner 666 667 .seealso: TSGLSetType(), TSGLRegister() 668 J*/ 669 typedef const char* TSGLType; 670 #define TSGL_IRKS "irks" 671 672 /*MC 673 TSGLRegisterDynamic - adds a TSGL implementation 674 675 Synopsis: 676 PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 677 678 Not Collective 679 680 Input Parameters: 681 + name_scheme - name of user-defined general linear scheme 682 . path - path (either absolute or relative) the library containing this scheme 683 . name_create - name of routine to create method context 684 - routine_create - routine to create method context 685 686 Notes: 687 TSGLRegisterDynamic() may be called multiple times to add several user-defined families. 688 689 If dynamic libraries are used, then the fourth input argument (routine_create) 690 is ignored. 691 692 Sample usage: 693 .vb 694 TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 695 "MySchemeCreate",MySchemeCreate); 696 .ve 697 698 Then, your scheme can be chosen with the procedural interface via 699 $ TSGLSetType(ts,"my_scheme") 700 or at runtime via the option 701 $ -ts_gl_type my_scheme 702 703 Level: advanced 704 705 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 706 and others of the form ${any_environmental_variable} occuring in pathname will be 707 replaced with appropriate values. 708 709 .keywords: TSGL, register 710 711 .seealso: TSGLRegisterAll() 712 M*/ 713 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 714 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 715 #else 716 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 717 #endif 718 719 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 720 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(const char[]); 721 PETSC_EXTERN PetscErrorCode TSGLRegisterDestroy(void); 722 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(const char[]); 723 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void); 724 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType); 725 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 726 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType); 727 728 /*J 729 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 730 731 Level: beginner 732 733 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 734 J*/ 735 typedef const char* TSARKIMEXType; 736 #define TSARKIMEX1BEE "1bee" 737 #define TSARKIMEXA2 "a2" 738 #define TSARKIMEXL2 "l2" 739 #define TSARKIMEXARS122 "ars122" 740 #define TSARKIMEX2C "2c" 741 #define TSARKIMEX2D "2d" 742 #define TSARKIMEX2E "2e" 743 #define TSARKIMEXPRSSP2 "prssp2" 744 #define TSARKIMEX3 "3" 745 #define TSARKIMEXBPR3 "bpr3" 746 #define TSARKIMEXARS443 "ars443" 747 #define TSARKIMEX4 "4" 748 #define TSARKIMEX5 "5" 749 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 750 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 751 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 752 PETSC_EXTERN PetscErrorCode TSARKIMEXRegister(TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]); 753 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 754 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(const char path[]); 755 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 756 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void); 757 758 /*J 759 TSRosWType - String with the name of a Rosenbrock-W method. 760 761 Level: beginner 762 763 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 764 J*/ 765 typedef const char* TSRosWType; 766 #define TSROSW2M "2m" 767 #define TSROSW2P "2p" 768 #define TSROSWRA3PW "ra3pw" 769 #define TSROSWRA34PW2 "ra34pw2" 770 #define TSROSWRODAS3 "rodas3" 771 #define TSROSWSANDU3 "sandu3" 772 #define TSROSWASSP3P3S1C "assp3p3s1c" 773 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 774 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 775 #define TSROSWARK3 "ark3" 776 #define TSROSWTHETA1 "theta1" 777 #define TSROSWTHETA2 "theta2" 778 #define TSROSWGRK4T "grk4t" 779 #define TSROSWSHAMP4 "shamp4" 780 #define TSROSWVELDD4 "veldd4" 781 #define TSROSW4L "4l" 782 783 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*); 784 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType); 785 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 786 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 787 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 788 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 789 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(const char path[]); 790 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 791 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void); 792 793 /* 794 PETSc interface to Sundials 795 */ 796 #ifdef PETSC_HAVE_SUNDIALS 797 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 798 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 799 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 800 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 801 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 802 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 803 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 804 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 805 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 806 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 807 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 808 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 809 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 810 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 811 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 812 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 813 #endif 814 815 PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal); 816 817 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 818 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 819 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 820 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 821 822 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 823 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 824 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 825 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 826 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 827 828 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 829 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 830 831 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 832 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 833 834 #endif 835