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, TSDestroy() 17 S*/ 18 typedef struct _p_TS* TS; 19 20 /*J 21 TSType - String with the name of a PETSc TS method. 22 23 Level: beginner 24 25 .seealso: TSSetType(), TS, TSRegister() 26 J*/ 27 typedef const char* TSType; 28 #define TSEULER "euler" 29 #define TSBEULER "beuler" 30 #define TSPSEUDO "pseudo" 31 #define TSCN "cn" 32 #define TSSUNDIALS "sundials" 33 #define TSRK "rk" 34 #define TSPYTHON "python" 35 #define TSTHETA "theta" 36 #define TSALPHA "alpha" 37 #define TSGL "gl" 38 #define TSSSP "ssp" 39 #define TSARKIMEX "arkimex" 40 #define TSROSW "rosw" 41 #define TSEIMEX "eimex" 42 #define TSMIMEX "mimex" 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 /*E 53 TSEquationType - type of TS problem that is solved 54 55 Level: beginner 56 57 Developer Notes: this must match petsc/finclude/petscts.h 58 59 Supported types are: 60 TS_EQ_UNSPECIFIED (default) 61 TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0 62 TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0 63 64 .seealso: TSGetEquationType(), TSSetEquationType() 65 E*/ 66 typedef enum { 67 TS_EQ_UNSPECIFIED = -1, 68 TS_EQ_EXPLICIT = 0, 69 TS_EQ_ODE_EXPLICIT = 1, 70 TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100, 71 TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200, 72 TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300, 73 TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500, 74 TS_EQ_IMPLICIT = 1000, 75 TS_EQ_ODE_IMPLICIT = 1001, 76 TS_EQ_DAE_IMPLICIT_INDEX1 = 1100, 77 TS_EQ_DAE_IMPLICIT_INDEX2 = 1200, 78 TS_EQ_DAE_IMPLICIT_INDEX3 = 1300, 79 TS_EQ_DAE_IMPLICIT_INDEXHI = 1500 80 } TSEquationType; 81 PETSC_EXTERN const char *const*TSEquationTypes; 82 83 /*E 84 TSConvergedReason - reason a TS method has converged or not 85 86 Level: beginner 87 88 Developer Notes: this must match petsc/finclude/petscts.h 89 90 Each reason has its own manual page. 91 92 .seealso: TSGetConvergedReason() 93 E*/ 94 typedef enum { 95 TS_CONVERGED_ITERATING = 0, 96 TS_CONVERGED_TIME = 1, 97 TS_CONVERGED_ITS = 2, 98 TS_CONVERGED_USER = 3, 99 TS_CONVERGED_EVENT = 4, 100 TS_DIVERGED_NONLINEAR_SOLVE = -1, 101 TS_DIVERGED_STEP_REJECTED = -2 102 } TSConvergedReason; 103 PETSC_EXTERN const char *const*TSConvergedReasons; 104 105 /*MC 106 TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve() 107 108 Level: beginner 109 110 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt() 111 M*/ 112 113 /*MC 114 TS_CONVERGED_TIME - the final time was reached 115 116 Level: beginner 117 118 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime() 119 M*/ 120 121 /*MC 122 TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time 123 124 Level: beginner 125 126 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration() 127 M*/ 128 129 /*MC 130 TS_CONVERGED_USER - user requested termination 131 132 Level: beginner 133 134 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 135 M*/ 136 137 /*MC 138 TS_CONVERGED_EVENT - user requested termination on event detection 139 140 Level: beginner 141 142 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 143 M*/ 144 145 /*MC 146 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 147 148 Level: beginner 149 150 Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures. 151 152 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures() 153 M*/ 154 155 /*MC 156 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 157 158 Level: beginner 159 160 Notes: See TSSetMaxStepRejections() for how to allow more step rejections. 161 162 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections() 163 M*/ 164 165 /*E 166 TSExactFinalTimeOption - option for handling of final time step 167 168 Level: beginner 169 170 Developer Notes: this must match petsc/finclude/petscts.h 171 172 $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded 173 $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time 174 $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time 175 .seealso: TSGetConvergedReason(), TSSetExactFinalTime() 176 177 E*/ 178 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption; 179 PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 180 181 182 /* Logging support */ 183 PETSC_EXTERN PetscClassId TS_CLASSID; 184 PETSC_EXTERN PetscClassId DMTS_CLASSID; 185 186 PETSC_EXTERN PetscErrorCode TSInitializePackage(void); 187 188 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*); 189 PETSC_EXTERN PetscErrorCode TSClone(TS,TS*); 190 PETSC_EXTERN PetscErrorCode TSDestroy(TS*); 191 192 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType); 193 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*); 194 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 195 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 196 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 197 198 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 199 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 200 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 201 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 202 PETSC_EXTERN PetscErrorCode TSSetUp(TS); 203 PETSC_EXTERN PetscErrorCode TSReset(TS); 204 205 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec); 206 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*); 207 208 /*S 209 TSTrajectory - Abstract PETSc object that storing the trajectory (solution of ODE/ADE at each time step and stage) 210 211 Level: advanced 212 213 Concepts: ODE solvers, adjoint 214 215 .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy() 216 S*/ 217 typedef struct _p_TSTrajectory* TSTrajectory; 218 219 /*J 220 TSTrajectoryType - String with the name of a PETSc TS trajectory storage method 221 222 Level: intermediate 223 224 .seealso: TSSetType(), TS, TSRegister(), TSTrajectoryCreate(), TSTrajectorySetType() 225 J*/ 226 typedef const char* TSTrajectoryType; 227 #define TSTRAJECTORYBASIC "basic" 228 #define TSTRAJECTORYSINGLEFILE "singlefile" 229 #define TSTRAJECTORYMEMORY "memory" 230 231 PETSC_EXTERN PetscFunctionList TSTrajectoryList; 232 PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID; 233 PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled; 234 235 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS); 236 237 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*); 238 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*); 239 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,const TSTrajectoryType); 240 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec); 241 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*); 242 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS); 243 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void); 244 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS); 245 246 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*); 247 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**); 248 PETSC_EXTERN PetscErrorCode TSSetCostIntegrand(TS,PetscInt,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),void*); 249 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*); 250 251 PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*); 252 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS); 253 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt); 254 255 PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat); 256 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS); 257 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS); 258 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*); 259 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*); 260 PETSC_EXTERN PetscErrorCode TSAdjointComputeCostIntegrand(TS,PetscReal,Vec,Vec); 261 262 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 263 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 264 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption); 265 266 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 267 268 typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx; 269 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *); 270 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*); 271 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*); 272 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*); 273 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*); 274 275 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*); 276 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 277 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 278 279 PETSC_EXTERN PetscErrorCode TSStep(TS); 280 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 281 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec); 282 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*); 283 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType); 284 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 285 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason); 286 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*); 287 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*); 288 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*); 289 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*); 290 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt); 291 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*); 292 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt); 293 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool); 294 PETSC_EXTERN PetscErrorCode TSRollBack(TS); 295 PETSC_EXTERN PetscErrorCode TSGetTotalSteps(TS,PetscInt*); 296 297 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**); 298 299 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 300 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*); 301 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*); 302 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal); 303 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 304 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal); 305 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*); 306 307 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 308 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*); 309 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 310 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 311 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 312 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 313 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool); 314 315 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*); 316 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*); 317 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*); 318 319 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 320 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 321 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 322 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 323 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 324 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 325 326 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 327 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*); 328 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 329 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 330 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec); 331 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec); 332 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 333 334 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 335 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal)); 336 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*)); 337 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 338 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 339 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal); 340 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*); 341 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 342 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool); 343 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 344 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 345 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*); 346 PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*); 347 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*); 348 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*); 349 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 350 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*); 351 352 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 353 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*); 354 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 355 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 356 357 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 358 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *); 359 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 360 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 361 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 362 363 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]); 364 365 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 366 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat); 367 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 368 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool); 369 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*); 370 371 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 372 373 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*); 374 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**); 375 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*); 376 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**); 377 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*); 378 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**); 379 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*); 380 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**); 381 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*); 382 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**); 383 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*); 384 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**); 385 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*); 386 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal); 387 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, const PetscReal[], PetscInt, PetscScalar *, void *), void **); 388 389 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*); 390 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*); 391 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*); 392 393 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 394 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 395 396 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 397 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*); 398 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 399 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*); 400 401 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *); 402 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *); 403 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *); 404 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *); 405 406 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *); 407 408 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 409 typedef struct { 410 Vec ray; 411 VecScatter scatter; 412 PetscViewer viewer; 413 TSMonitorLGCtx lgctx; 414 } TSMonitorDMDARayCtx; 415 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**); 416 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*); 417 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*); 418 419 420 /* Dynamic creation and loading functions */ 421 PETSC_EXTERN PetscFunctionList TSList; 422 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*); 423 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType); 424 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS)); 425 426 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 427 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES); 428 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 429 430 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 431 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer); 432 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 433 434 435 #define TS_FILE_CLASSID 1211225 436 437 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 438 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 439 440 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 441 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 442 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 443 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 444 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*); 445 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **); 446 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *); 447 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*); 448 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*); 449 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*); 450 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *); 451 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 452 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 453 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 454 455 typedef struct _n_TSMonitorEnvelopeCtx* TSMonitorEnvelopeCtx; 456 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*); 457 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*); 458 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*); 459 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*); 460 461 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 462 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 463 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 464 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 465 466 PETSC_EXTERN PetscErrorCode TSSetEventMonitor(TS,PetscInt,PetscInt*,PetscBool*,PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*); 467 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal*); 468 /*J 469 TSSSPType - string with the name of TSSSP scheme. 470 471 Level: beginner 472 473 .seealso: TSSSPSetType(), TS 474 J*/ 475 typedef const char* TSSSPType; 476 #define TSSSPRKS2 "rks2" 477 #define TSSSPRKS3 "rks3" 478 #define TSSSPRK104 "rk104" 479 480 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 481 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 482 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 483 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 484 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void); 485 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void); 486 487 /*S 488 TSAdapt - Abstract object that manages time-step adaptivity 489 490 Level: beginner 491 492 .seealso: TS, TSAdaptCreate(), TSAdaptType 493 S*/ 494 typedef struct _p_TSAdapt *TSAdapt; 495 496 /*E 497 TSAdaptType - String with the name of TSAdapt scheme. 498 499 Level: beginner 500 501 .seealso: TSAdaptSetType(), TS 502 E*/ 503 typedef const char *TSAdaptType; 504 #define TSADAPTBASIC "basic" 505 #define TSADAPTNONE "none" 506 #define TSADAPTCFL "cfl" 507 508 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*); 509 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt)); 510 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void); 511 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 512 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 513 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 514 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 515 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 516 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 517 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 518 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 519 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*); 520 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 521 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer); 522 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptions*,TSAdapt); 523 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt); 524 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 525 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 526 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 527 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*)); 528 529 /*S 530 TSGLAdapt - Abstract object that manages time-step adaptivity 531 532 Level: beginner 533 534 Developer Notes: 535 This functionality should be replaced by the TSAdapt. 536 537 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 538 S*/ 539 typedef struct _p_TSGLAdapt *TSGLAdapt; 540 541 /*J 542 TSGLAdaptType - String with the name of TSGLAdapt scheme 543 544 Level: beginner 545 546 .seealso: TSGLAdaptSetType(), TS 547 J*/ 548 typedef const char *TSGLAdaptType; 549 #define TSGLADAPT_NONE "none" 550 #define TSGLADAPT_SIZE "size" 551 #define TSGLADAPT_BOTH "both" 552 553 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt)); 554 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void); 555 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void); 556 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 557 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType); 558 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 559 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 560 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 561 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptions*,TSGLAdapt); 562 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 563 564 /*J 565 TSGLAcceptType - String with the name of TSGLAccept scheme 566 567 Level: beginner 568 569 .seealso: TSGLSetAcceptType(), TS 570 J*/ 571 typedef const char *TSGLAcceptType; 572 #define TSGLACCEPT_ALWAYS "always" 573 574 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 575 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction); 576 577 /*J 578 TSGLType - family of time integration method within the General Linear class 579 580 Level: beginner 581 582 .seealso: TSGLSetType(), TSGLRegister() 583 J*/ 584 typedef const char* TSGLType; 585 #define TSGL_IRKS "irks" 586 587 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS)); 588 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void); 589 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void); 590 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType); 591 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 592 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType); 593 594 /*J 595 TSEIMEXType - String with the name of an Extrapolated IMEX method. 596 597 Level: beginner 598 599 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister() 600 J*/ 601 #define TSEIMEXType char* 602 603 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt); 604 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt); 605 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool); 606 607 /*J 608 TSRKType - String with the name of a Runge-Kutta method. 609 610 Level: beginner 611 612 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister() 613 J*/ 614 typedef const char* TSRKType; 615 #define TSRK1FE "1fe" 616 #define TSRK2A "2a" 617 #define TSRK3 "3" 618 #define TSRK3BS "3bs" 619 #define TSRK4 "4" 620 #define TSRK5F "5f" 621 #define TSRK5DP "5dp" 622 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*); 623 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType); 624 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool); 625 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 626 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void); 627 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void); 628 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void); 629 630 /*J 631 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 632 633 Level: beginner 634 635 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 636 J*/ 637 typedef const char* TSARKIMEXType; 638 #define TSARKIMEX1BEE "1bee" 639 #define TSARKIMEXA2 "a2" 640 #define TSARKIMEXL2 "l2" 641 #define TSARKIMEXARS122 "ars122" 642 #define TSARKIMEX2C "2c" 643 #define TSARKIMEX2D "2d" 644 #define TSARKIMEX2E "2e" 645 #define TSARKIMEXPRSSP2 "prssp2" 646 #define TSARKIMEX3 "3" 647 #define TSARKIMEXBPR3 "bpr3" 648 #define TSARKIMEXARS443 "ars443" 649 #define TSARKIMEX4 "4" 650 #define TSARKIMEX5 "5" 651 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 652 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 653 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 654 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[]); 655 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 656 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void); 657 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 658 659 /*J 660 TSRosWType - String with the name of a Rosenbrock-W method. 661 662 Level: beginner 663 664 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 665 J*/ 666 typedef const char* TSRosWType; 667 #define TSROSW2M "2m" 668 #define TSROSW2P "2p" 669 #define TSROSWRA3PW "ra3pw" 670 #define TSROSWRA34PW2 "ra34pw2" 671 #define TSROSWRODAS3 "rodas3" 672 #define TSROSWSANDU3 "sandu3" 673 #define TSROSWASSP3P3S1C "assp3p3s1c" 674 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 675 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 676 #define TSROSWARK3 "ark3" 677 #define TSROSWTHETA1 "theta1" 678 #define TSROSWTHETA2 "theta2" 679 #define TSROSWGRK4T "grk4t" 680 #define TSROSWSHAMP4 "shamp4" 681 #define TSROSWVELDD4 "veldd4" 682 #define TSROSW4L "4l" 683 684 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*); 685 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType); 686 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 687 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 688 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 689 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 690 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void); 691 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 692 693 /* 694 PETSc interface to Sundials 695 */ 696 #ifdef PETSC_HAVE_SUNDIALS 697 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 698 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 699 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 700 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 701 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 702 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 703 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 704 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 705 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 706 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 707 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 708 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 709 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 710 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 711 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 712 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 713 #endif 714 715 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 716 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 717 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 718 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 719 720 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 721 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 722 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 723 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 724 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 725 726 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 727 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 728 729 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 730 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*); 731 732 #endif 733