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 PetscFunctionList 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 #include "petscts.h" 335 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 336 337 Not Collective 338 339 Input Parameters: 340 + name - The name of a new user-defined creation routine 341 . path - The path (either absolute or relative) of the library containing this routine 342 . func_name - The name of the creation routine 343 - create_func - The creation routine itself 344 345 Notes: 346 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 347 348 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 349 350 Sample usage: 351 .vb 352 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 353 .ve 354 355 Then, your ts type can be chosen with the procedural interface via 356 .vb 357 TS ts; 358 TSCreate(MPI_Comm, &ts); 359 TSSetType(ts, "my_ts") 360 .ve 361 or at runtime via the option 362 .vb 363 -ts_type my_ts 364 .ve 365 366 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 367 If your function is not being put into a shared library then use TSRegister() instead 368 369 Level: advanced 370 371 .keywords: TS, register 372 .seealso: TSRegisterAll(), TSRegisterDestroy() 373 M*/ 374 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 375 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 376 #else 377 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 378 #endif 379 380 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 381 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES); 382 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 383 384 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 385 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer); 386 387 #define TS_FILE_CLASSID 1211225 388 389 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 390 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 391 392 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 393 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 394 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 395 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 396 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 397 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 398 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 399 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 400 401 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 402 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 403 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 404 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 405 406 /*J 407 TSSSPType - string with the name of TSSSP scheme. 408 409 Level: beginner 410 411 .seealso: TSSSPSetType(), TS 412 J*/ 413 typedef const char* TSSSPType; 414 #define TSSSPRKS2 "rks2" 415 #define TSSSPRKS3 "rks3" 416 #define TSSSPRK104 "rk104" 417 418 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 419 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 420 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 421 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 422 423 /*S 424 TSAdapt - Abstract object that manages time-step adaptivity 425 426 Level: beginner 427 428 .seealso: TS, TSAdaptCreate(), TSAdaptType 429 S*/ 430 typedef struct _p_TSAdapt *TSAdapt; 431 432 /*E 433 TSAdaptType - String with the name of TSAdapt scheme or the creation function 434 with an optional dynamic library name, for example 435 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 436 437 Level: beginner 438 439 .seealso: TSAdaptSetType(), TS 440 E*/ 441 typedef const char *TSAdaptType; 442 #define TSADAPTBASIC "basic" 443 #define TSADAPTNONE "none" 444 #define TSADAPTCFL "cfl" 445 446 /*MC 447 TSAdaptRegisterDynamic - adds a TSAdapt implementation 448 449 Synopsis: 450 #include "petscts.h" 451 PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 452 453 Not Collective 454 455 Input Parameters: 456 + name_scheme - name of user-defined adaptivity scheme 457 . path - path (either absolute or relative) the library containing this scheme 458 . name_create - name of routine to create method context 459 - routine_create - routine to create method context 460 461 Notes: 462 TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 463 464 If dynamic libraries are used, then the fourth input argument (routine_create) 465 is ignored. 466 467 Sample usage: 468 .vb 469 TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 470 "MySchemeCreate",MySchemeCreate); 471 .ve 472 473 Then, your scheme can be chosen with the procedural interface via 474 $ TSAdaptSetType(ts,"my_scheme") 475 or at runtime via the option 476 $ -ts_adapt_type my_scheme 477 478 Level: advanced 479 480 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 481 and others of the form ${any_environmental_variable} occuring in pathname will be 482 replaced with appropriate values. 483 484 .keywords: TSAdapt, register 485 486 .seealso: TSAdaptRegisterAll() 487 M*/ 488 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 489 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,0) 490 #else 491 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,d) 492 #endif 493 494 PETSC_EXTERN PetscErrorCode TSGetTSAdapt(TS,TSAdapt*); 495 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt)); 496 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(const char[]); 497 PETSC_EXTERN PetscErrorCode TSAdaptRegisterDestroy(void); 498 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(const char[]); 499 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 500 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 501 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 502 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 503 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 504 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 505 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 506 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 507 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*); 508 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 509 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer); 510 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt); 511 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 512 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 513 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 514 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*)); 515 516 /*S 517 TSGLAdapt - Abstract object that manages time-step adaptivity 518 519 Level: beginner 520 521 Developer Notes: 522 This functionality should be replaced by the TSAdapt. 523 524 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 525 S*/ 526 typedef struct _p_TSGLAdapt *TSGLAdapt; 527 528 /*J 529 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 530 with an optional dynamic library name, for example 531 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 532 533 Level: beginner 534 535 .seealso: TSGLAdaptSetType(), TS 536 J*/ 537 typedef const char *TSGLAdaptType; 538 #define TSGLADAPT_NONE "none" 539 #define TSGLADAPT_SIZE "size" 540 #define TSGLADAPT_BOTH "both" 541 542 /*MC 543 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 544 545 Synopsis: 546 #include "petscts.h" 547 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 548 549 Not Collective 550 551 Input Parameters: 552 + name_scheme - name of user-defined adaptivity scheme 553 . path - path (either absolute or relative) the library containing this scheme 554 . name_create - name of routine to create method context 555 - routine_create - routine to create method context 556 557 Notes: 558 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 559 560 If dynamic libraries are used, then the fourth input argument (routine_create) 561 is ignored. 562 563 Sample usage: 564 .vb 565 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 566 "MySchemeCreate",MySchemeCreate); 567 .ve 568 569 Then, your scheme can be chosen with the procedural interface via 570 $ TSGLAdaptSetType(ts,"my_scheme") 571 or at runtime via the option 572 $ -ts_adapt_type my_scheme 573 574 Level: advanced 575 576 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 577 and others of the form ${any_environmental_variable} occuring in pathname will be 578 replaced with appropriate values. 579 580 .keywords: TSGLAdapt, register 581 582 .seealso: TSGLAdaptRegisterAll() 583 M*/ 584 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 585 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 586 #else 587 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 588 #endif 589 590 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 591 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(const char[]); 592 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterDestroy(void); 593 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(const char[]); 594 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void); 595 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 596 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType); 597 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 598 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 599 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 600 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 601 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 602 603 /*J 604 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 605 with an optional dynamic library name, for example 606 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 607 608 Level: beginner 609 610 .seealso: TSGLSetAcceptType(), TS 611 J*/ 612 typedef const char *TSGLAcceptType; 613 #define TSGLACCEPT_ALWAYS "always" 614 615 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 616 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 617 618 /*MC 619 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 620 621 Synopsis: 622 #include "petscts.h" 623 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 624 625 Not Collective 626 627 Input Parameters: 628 + name_scheme - name of user-defined acceptance scheme 629 . path - path (either absolute or relative) the library containing this scheme 630 . name_create - name of routine to create method context 631 - routine_create - routine to create method context 632 633 Notes: 634 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 635 636 If dynamic libraries are used, then the fourth input argument (routine_create) 637 is ignored. 638 639 Sample usage: 640 .vb 641 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 642 "MySchemeCreate",MySchemeCreate); 643 .ve 644 645 Then, your scheme can be chosen with the procedural interface via 646 $ TSGLSetAcceptType(ts,"my_scheme") 647 or at runtime via the option 648 $ -ts_gl_accept_type my_scheme 649 650 Level: advanced 651 652 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 653 and others of the form ${any_environmental_variable} occuring in pathname will be 654 replaced with appropriate values. 655 656 .keywords: TSGL, TSGLAcceptType, register 657 658 .seealso: TSGLRegisterAll() 659 M*/ 660 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 661 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 662 #else 663 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 664 #endif 665 666 /*J 667 TSGLType - family of time integration method within the General Linear class 668 669 Level: beginner 670 671 .seealso: TSGLSetType(), TSGLRegister() 672 J*/ 673 typedef const char* TSGLType; 674 #define TSGL_IRKS "irks" 675 676 /*MC 677 TSGLRegisterDynamic - adds a TSGL implementation 678 679 Synopsis: 680 #include "petscts.h" 681 PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 682 683 Not Collective 684 685 Input Parameters: 686 + name_scheme - name of user-defined general linear scheme 687 . path - path (either absolute or relative) the library containing this scheme 688 . name_create - name of routine to create method context 689 - routine_create - routine to create method context 690 691 Notes: 692 TSGLRegisterDynamic() may be called multiple times to add several user-defined families. 693 694 If dynamic libraries are used, then the fourth input argument (routine_create) 695 is ignored. 696 697 Sample usage: 698 .vb 699 TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 700 "MySchemeCreate",MySchemeCreate); 701 .ve 702 703 Then, your scheme can be chosen with the procedural interface via 704 $ TSGLSetType(ts,"my_scheme") 705 or at runtime via the option 706 $ -ts_gl_type my_scheme 707 708 Level: advanced 709 710 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 711 and others of the form ${any_environmental_variable} occuring in pathname will be 712 replaced with appropriate values. 713 714 .keywords: TSGL, register 715 716 .seealso: TSGLRegisterAll() 717 M*/ 718 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 719 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 720 #else 721 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 722 #endif 723 724 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 725 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(const char[]); 726 PETSC_EXTERN PetscErrorCode TSGLRegisterDestroy(void); 727 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(const char[]); 728 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void); 729 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType); 730 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 731 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType); 732 733 /*J 734 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 735 736 Level: beginner 737 738 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 739 J*/ 740 typedef const char* TSARKIMEXType; 741 #define TSARKIMEX1BEE "1bee" 742 #define TSARKIMEXA2 "a2" 743 #define TSARKIMEXL2 "l2" 744 #define TSARKIMEXARS122 "ars122" 745 #define TSARKIMEX2C "2c" 746 #define TSARKIMEX2D "2d" 747 #define TSARKIMEX2E "2e" 748 #define TSARKIMEXPRSSP2 "prssp2" 749 #define TSARKIMEX3 "3" 750 #define TSARKIMEXBPR3 "bpr3" 751 #define TSARKIMEXARS443 "ars443" 752 #define TSARKIMEX4 "4" 753 #define TSARKIMEX5 "5" 754 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 755 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 756 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 757 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[]); 758 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 759 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(const char path[]); 760 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 761 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void); 762 763 /*J 764 TSRosWType - String with the name of a Rosenbrock-W method. 765 766 Level: beginner 767 768 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 769 J*/ 770 typedef const char* TSRosWType; 771 #define TSROSW2M "2m" 772 #define TSROSW2P "2p" 773 #define TSROSWRA3PW "ra3pw" 774 #define TSROSWRA34PW2 "ra34pw2" 775 #define TSROSWRODAS3 "rodas3" 776 #define TSROSWSANDU3 "sandu3" 777 #define TSROSWASSP3P3S1C "assp3p3s1c" 778 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 779 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 780 #define TSROSWARK3 "ark3" 781 #define TSROSWTHETA1 "theta1" 782 #define TSROSWTHETA2 "theta2" 783 #define TSROSWGRK4T "grk4t" 784 #define TSROSWSHAMP4 "shamp4" 785 #define TSROSWVELDD4 "veldd4" 786 #define TSROSW4L "4l" 787 788 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*); 789 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType); 790 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 791 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 792 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 793 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 794 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(const char path[]); 795 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 796 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void); 797 798 /* 799 PETSc interface to Sundials 800 */ 801 #ifdef PETSC_HAVE_SUNDIALS 802 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 803 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 804 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 805 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 806 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 807 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 808 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 809 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 810 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 811 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 812 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 813 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 814 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 815 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 816 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 817 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 818 #endif 819 820 PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal); 821 822 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 823 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 824 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 825 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 826 827 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 828 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 829 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 830 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 831 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 832 833 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 834 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 835 836 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 837 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 838 839 #endif 840