Actual source code: petscts.h

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: /*
  2:    User interface for the timestepping package. This package
  3:    is for use in solving time-dependent PDEs.
  4: */
  7: #include <petscsnes.h>

  9: /*S
 10:      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)

 12:    Level: beginner

 14:   Concepts: ODE solvers

 16: .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
 17: S*/
 18: typedef struct _p_TS* TS;

 20: /*J
 21:     TSType - String with the name of a PETSc TS method.

 23:    Level: beginner

 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

 46:    Level: beginner

 48: .seealso: TSCreate()
 49: E*/
 50: typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;

 52: /*E
 53:    TSEquationType - type of TS problem that is solved

 55:    Level: beginner

 57:    Developer Notes: this must match petsc/finclude/petscts.h

 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

 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;

 83: /*E
 84:    TSConvergedReason - reason a TS method has converged or not

 86:    Level: beginner

 88:    Developer Notes: this must match petsc/finclude/petscts.h

 90:    Each reason has its own manual page.

 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;

105: /*MC
106:    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()

108:    Level: beginner

110: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
111: M*/

113: /*MC
114:    TS_CONVERGED_TIME - the final time was reached

116:    Level: beginner

118: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime()
119: M*/

121: /*MC
122:    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time

124:    Level: beginner

126: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration()
127: M*/

129: /*MC
130:    TS_CONVERGED_USER - user requested termination

132:    Level: beginner

134: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
135: M*/

137: /*MC
138:    TS_CONVERGED_EVENT - user requested termination on event detection

140:    Level: beginner

142: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
143: M*/

145: /*MC
146:    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed

148:    Level: beginner

150:    Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.

152: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
153: M*/

155: /*MC
156:    TS_DIVERGED_STEP_REJECTED - too many steps were rejected

158:    Level: beginner

160:    Notes: See TSSetMaxStepRejections() for how to allow more step rejections.

162: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
163: M*/

165: /*E
166:    TSExactFinalTimeOption - option for handling of final time step

168:    Level: beginner

170:    Developer Notes: this must match petsc/finclude/petscts.h

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()

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[];


182: /* Logging support */
183: PETSC_EXTERN PetscClassId TS_CLASSID;
184: PETSC_EXTERN PetscClassId DMTS_CLASSID;

186: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);

188: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
189: PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
190: PETSC_EXTERN PetscErrorCode TSDestroy(TS*);

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);

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);

205: PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
206: PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);

208: /*S
209:      TSTrajectory - Abstract PETSc object that storing the trajectory (solution of ODE/ADE at each time step and stage)

211:    Level: advanced

213:   Concepts: ODE solvers, adjoint

215: .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
216: S*/
217: typedef struct _p_TSTrajectory* TSTrajectory;

219: /*J
220:     TSTrajectoryType - String with the name of a PETSc TS trajectory storage method

222:    Level: intermediate

224: .seealso: TSSetType(), TS, TSRegister(), TSTrajectoryCreate(), TSTrajectorySetType()
225: J*/
226: typedef const char* TSTrajectoryType;
227: #define TSTRAJECTORYBASIC      "basic"
228: #define TSTRAJECTORYSINGLEFILE "singlefile"

230: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
231: PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
232: PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;

234: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);

236: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
237: PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
238: PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,const TSTrajectoryType);
239: PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
240: PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
241: PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory);
242: PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);

244: PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
245: PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
246: 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*);
247: PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);

249: PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
250: PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
251: PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);

253: PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
254: PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
255: PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
256: PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
257: PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
258: PETSC_EXTERN PetscErrorCode TSAdjointComputeCostIntegrand(TS,PetscReal,Vec,Vec);

260: PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
261: PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
262: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);

264: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);

266: typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
267: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
268: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
269: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
270: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
271: PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);

273: PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
274: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
275: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);

277: PETSC_EXTERN PetscErrorCode TSStep(TS);
278: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
279: PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
280: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
281: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
282: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
283: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
284: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
285: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
286: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
287: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
288: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
289: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
290: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
291: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
292: PETSC_EXTERN PetscErrorCode TSRollBack(TS);
293: PETSC_EXTERN PetscErrorCode TSGetTotalSteps(TS,PetscInt*);

295: PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);

297: PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
298: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
299: PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
300: PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
301: PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
302: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
303: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);

305: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
306: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
307: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
308: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
309: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
310: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
311: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);

313: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
314: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
315: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);

317: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
318: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
319: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
320: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
321: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
322: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);

324: PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
325: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
326: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
327: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
328: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
329: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
330: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);

332: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
333: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
334: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
335: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
336: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
337: PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
338: PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
339: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
340: PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
341: PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
342: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
343: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
344: PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*);
345: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*);
346: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*);
347: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
348: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);

350: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
351: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
352: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
353: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);

355: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
356: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
357: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
358: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
359: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);

361: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);

363: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
364: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
365: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
366: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
367: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);

369: PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);

371: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
372: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
373: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
374: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
375: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
376: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
377: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
378: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
379: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
380: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
381: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
382: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**);
383: PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
384: PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
385: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, const PetscReal[], PetscInt, PetscScalar *, void *), void **);

387: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
388: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
389: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);

391: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
392: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));

394: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
395: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
396: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
397: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);

399: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
400: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
401: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
402: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);

404: PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);

406: typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
407: typedef struct {
408:   Vec            ray;
409:   VecScatter     scatter;
410:   PetscViewer    viewer;
411:   TSMonitorLGCtx lgctx;
412: } TSMonitorDMDARayCtx;
413: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
414: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
415: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);


418: /* Dynamic creation and loading functions */
419: PETSC_EXTERN PetscFunctionList TSList;
420: PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
421: PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
422: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));

424: PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
425: PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
426: PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);

428: PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
429: PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
430: PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}


433: #define TS_FILE_CLASSID 1211225

435: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
436: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);

438: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
439: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
440: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
441: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
442: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
443: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
444: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
445: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
446: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
447: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
448: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
449: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
450: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
451: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);

453: typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
454: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
455: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
456: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
457: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);

459: typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
460: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
461: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
462: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);

464: PETSC_EXTERN PetscErrorCode TSSetEventMonitor(TS,PetscInt,PetscInt*,PetscBool*,PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
465: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal*);
466: /*J
467:    TSSSPType - string with the name of TSSSP scheme.

469:    Level: beginner

471: .seealso: TSSSPSetType(), TS
472: J*/
473: typedef const char* TSSSPType;
474: #define TSSSPRKS2  "rks2"
475: #define TSSSPRKS3  "rks3"
476: #define TSSSPRK104 "rk104"

478: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
479: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
480: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
481: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
482: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
483: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);

485: /*S
486:    TSAdapt - Abstract object that manages time-step adaptivity

488:    Level: beginner

490: .seealso: TS, TSAdaptCreate(), TSAdaptType
491: S*/
492: typedef struct _p_TSAdapt *TSAdapt;

494: /*E
495:     TSAdaptType - String with the name of TSAdapt scheme.

497:    Level: beginner

499: .seealso: TSAdaptSetType(), TS
500: E*/
501: typedef const char *TSAdaptType;
502: #define TSADAPTBASIC "basic"
503: #define TSADAPTNONE  "none"
504: #define TSADAPTCFL   "cfl"

506: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
507: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
508: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
509: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
510: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
511: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
512: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
513: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
514: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
515: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
516: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
517: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
518: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
519: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
520: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptions*,TSAdapt);
521: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
522: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
523: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
524: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
525: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*));

527: /*S
528:    TSGLAdapt - Abstract object that manages time-step adaptivity

530:    Level: beginner

532:    Developer Notes:
533:    This functionality should be replaced by the TSAdapt.

535: .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
536: S*/
537: typedef struct _p_TSGLAdapt *TSGLAdapt;

539: /*J
540:     TSGLAdaptType - String with the name of TSGLAdapt scheme

542:    Level: beginner

544: .seealso: TSGLAdaptSetType(), TS
545: J*/
546: typedef const char *TSGLAdaptType;
547: #define TSGLADAPT_NONE "none"
548: #define TSGLADAPT_SIZE "size"
549: #define TSGLADAPT_BOTH "both"

551: PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt));
552: PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void);
553: PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
554: PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
555: PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
556: PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
557: PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
558: PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
559: PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptions*,TSGLAdapt);
560: PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);

562: /*J
563:     TSGLAcceptType - String with the name of TSGLAccept scheme

565:    Level: beginner

567: .seealso: TSGLSetAcceptType(), TS
568: J*/
569: typedef const char *TSGLAcceptType;
570: #define TSGLACCEPT_ALWAYS "always"

572: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
573: PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction);

575: /*J
576:   TSGLType - family of time integration method within the General Linear class

578:   Level: beginner

580: .seealso: TSGLSetType(), TSGLRegister()
581: J*/
582: typedef const char* TSGLType;
583: #define TSGL_IRKS   "irks"

585: PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS));
586: PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void);
587: PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
588: PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
589: PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
590: PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);

592: /*J
593:     TSEIMEXType - String with the name of an Extrapolated IMEX method.

595:    Level: beginner

597: .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
598: J*/
599: #define TSEIMEXType   char*

601: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
602: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
603: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);

605: /*J
606:     TSRKType - String with the name of a Runge-Kutta method.

608:    Level: beginner

610: .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
611: J*/
612: typedef const char* TSRKType;
613: #define TSRK1FE   "1fe"
614: #define TSRK2A    "2a"
615: #define TSRK3     "3"
616: #define TSRK3BS   "3bs"
617: #define TSRK4     "4"
618: #define TSRK5F    "5f"
619: #define TSRK5DP   "5dp"
620: PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
621: PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
622: PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
623: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
624: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
625: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
626: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);

628: /*J
629:     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.

631:    Level: beginner

633: .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
634: J*/
635: typedef const char* TSARKIMEXType;
636: #define TSARKIMEX1BEE   "1bee"
637: #define TSARKIMEXA2     "a2"
638: #define TSARKIMEXL2     "l2"
639: #define TSARKIMEXARS122 "ars122"
640: #define TSARKIMEX2C     "2c"
641: #define TSARKIMEX2D     "2d"
642: #define TSARKIMEX2E     "2e"
643: #define TSARKIMEXPRSSP2 "prssp2"
644: #define TSARKIMEX3      "3"
645: #define TSARKIMEXBPR3   "bpr3"
646: #define TSARKIMEXARS443 "ars443"
647: #define TSARKIMEX4      "4"
648: #define TSARKIMEX5      "5"
649: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
650: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
651: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
652: 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[]);
653: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
654: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
655: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);

657: /*J
658:     TSRosWType - String with the name of a Rosenbrock-W method.

660:    Level: beginner

662: .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
663: J*/
664: typedef const char* TSRosWType;
665: #define TSROSW2M          "2m"
666: #define TSROSW2P          "2p"
667: #define TSROSWRA3PW       "ra3pw"
668: #define TSROSWRA34PW2     "ra34pw2"
669: #define TSROSWRODAS3      "rodas3"
670: #define TSROSWSANDU3      "sandu3"
671: #define TSROSWASSP3P3S1C  "assp3p3s1c"
672: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
673: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
674: #define TSROSWARK3        "ark3"
675: #define TSROSWTHETA1      "theta1"
676: #define TSROSWTHETA2      "theta2"
677: #define TSROSWGRK4T       "grk4t"
678: #define TSROSWSHAMP4      "shamp4"
679: #define TSROSWVELDD4      "veldd4"
680: #define TSROSW4L          "4l"

682: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
683: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
684: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
685: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
686: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
687: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
688: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
689: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);

691: /*
692:        PETSc interface to Sundials
693: */
694: #ifdef PETSC_HAVE_SUNDIALS
695: typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
696: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
697: typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
698: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
699: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
700: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
701: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
702: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
703: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
704: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
705: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
706: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
707: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
708: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
709: PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
710: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
711: #endif

713: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
714: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
715: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
716: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);

718: PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
719: PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
720: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
721: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
722: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);

724: PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
725: PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);

727: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
728: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);

730: #endif