Actual source code: petscts.h

petsc-3.11.4 2019-09-28
Report Typos and Errors
  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>

  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 TSBASICSYMPLECTIC "basicsymplectic"
 31: #define TSPSEUDO          "pseudo"
 32: #define TSCN              "cn"
 33: #define TSSUNDIALS        "sundials"
 34: #define TSRK              "rk"
 35: #define TSPYTHON          "python"
 36: #define TSTHETA           "theta"
 37: #define TSALPHA           "alpha"
 38: #define TSALPHA2          "alpha2"
 39: #define TSGLLE            "glle"
 40: #define TSGLEE            "glee"
 41: #define TSSSP             "ssp"
 42: #define TSARKIMEX         "arkimex"
 43: #define TSROSW            "rosw"
 44: #define TSEIMEX           "eimex"
 45: #define TSMIMEX           "mimex"
 46: #define TSBDF             "bdf"
 47: #define TSRADAU5          "radau5"

 49: /*E
 50:     TSProblemType - Determines the type of problem this TS object is to be used to solve

 52:    Level: beginner

 54: .seealso: TSCreate()
 55: E*/
 56: typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;

 58: /*E
 59:    TSEquationType - type of TS problem that is solved

 61:    Level: beginner

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

 66:    Supported types are:
 67:      TS_EQ_UNSPECIFIED (default)
 68:      TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
 69:      TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0

 71: .seealso: TSGetEquationType(), TSSetEquationType()
 72: E*/
 73: typedef enum {
 74:   TS_EQ_UNSPECIFIED               = -1,
 75:   TS_EQ_EXPLICIT                  = 0,
 76:   TS_EQ_ODE_EXPLICIT              = 1,
 77:   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
 78:   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
 79:   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
 80:   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
 81:   TS_EQ_IMPLICIT                  = 1000,
 82:   TS_EQ_ODE_IMPLICIT              = 1001,
 83:   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
 84:   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
 85:   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
 86:   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
 87: } TSEquationType;
 88: PETSC_EXTERN const char *const*TSEquationTypes;

 90: /*E
 91:    TSConvergedReason - reason a TS method has converged or not

 93:    Level: beginner

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

 98:    Each reason has its own manual page.

100: .seealso: TSGetConvergedReason()
101: E*/
102: typedef enum {
103:   TS_CONVERGED_ITERATING      = 0,
104:   TS_CONVERGED_TIME           = 1,
105:   TS_CONVERGED_ITS            = 2,
106:   TS_CONVERGED_USER           = 3,
107:   TS_CONVERGED_EVENT          = 4,
108:   TS_CONVERGED_PSEUDO_FATOL   = 5,
109:   TS_CONVERGED_PSEUDO_FRTOL   = 6,
110:   TS_DIVERGED_NONLINEAR_SOLVE = -1,
111:   TS_DIVERGED_STEP_REJECTED   = -2
112: } TSConvergedReason;
113: PETSC_EXTERN const char *const*TSConvergedReasons;

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

118:    Level: beginner

120: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
121: M*/

123: /*MC
124:    TS_CONVERGED_TIME - the final time was reached

126:    Level: beginner

128: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxTime(), TSGetMaxTime(), TSGetSolveTime()
129: M*/

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

134:    Level: beginner

136: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxSteps(), TSGetMaxSteps()
137: M*/

139: /*MC
140:    TS_CONVERGED_USER - user requested termination

142:    Level: beginner

144: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
145: M*/

147: /*MC
148:    TS_CONVERGED_EVENT - user requested termination on event detection

150:    Level: beginner

152: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
153: M*/

155: /*MC
156:    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for TSPSEUDO

158:    Level: beginner

160:    Options Database:
161: .   -ts_pseudo_frtol <rtol>

163: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FATOL
164: M*/

166: /*MC
167:    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for TSPSEUDO

169:    Level: beginner

171:    Options Database:
172: .   -ts_pseudo_fatol <atol>

174: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FRTOL
175: M*/

177: /*MC
178:    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed

180:    Level: beginner

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

185: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
186: M*/

188: /*MC
189:    TS_DIVERGED_STEP_REJECTED - too many steps were rejected

191:    Level: beginner

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

196: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
197: M*/

199: /*E
200:    TSExactFinalTimeOption - option for handling of final time step

202:    Level: beginner

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

207: $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
208: $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
209: $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time

211: .seealso: TSGetConvergedReason(), TSSetExactFinalTime(), TSGetExactFinalTime()

213: E*/
214: typedef enum {TS_EXACTFINALTIME_UNSPECIFIED=0,TS_EXACTFINALTIME_STEPOVER=1,TS_EXACTFINALTIME_INTERPOLATE=2,TS_EXACTFINALTIME_MATCHSTEP=3} TSExactFinalTimeOption;
215: PETSC_EXTERN const char *const TSExactFinalTimeOptions[];

217: /* Logging support */
218: PETSC_EXTERN PetscClassId TS_CLASSID;
219: PETSC_EXTERN PetscClassId DMTS_CLASSID;
220: PETSC_EXTERN PetscClassId TSADAPT_CLASSID;

222: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
223: PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);

225: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
226: PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
227: PETSC_EXTERN PetscErrorCode TSDestroy(TS*);

229: PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
230: PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
231: PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
232: PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
233: PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
234: PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);

236: PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
237: PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
238: PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
239: PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
240: PETSC_EXTERN PetscErrorCode TSSetUp(TS);
241: PETSC_EXTERN PetscErrorCode TSReset(TS);

243: PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
244: PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);

246: PETSC_EXTERN PetscErrorCode TS2SetSolution(TS,Vec,Vec);
247: PETSC_EXTERN PetscErrorCode TS2GetSolution(TS,Vec*,Vec*);

249: PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS,PetscInt*,Vec*);
250: PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS,Vec*);
251: PETSC_EXTERN PetscErrorCode TSGetTimeError(TS,PetscInt,Vec*);
252: PETSC_EXTERN PetscErrorCode TSSetTimeError(TS,Vec);

254: PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
255: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS,PetscReal,Vec,Mat);
256: PETSC_EXTERN PetscErrorCode TSComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
257: PETSC_EXTERN PetscErrorCode TSComputeDRDYFunction(TS,PetscReal,Vec,Vec*);

259: /*S
260:      TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)

262:    Level: advanced

264:   Concepts: ODE solvers, trajectory

266: .seealso:  TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectorySetType(), TSTrajectoryDestroy(), TSTrajectoryReset()
267: S*/
268: typedef struct _p_TSTrajectory* TSTrajectory;

270: /*J
271:     TSTrajectorySetType - String with the name of a PETSc TS trajectory storage method

273:    Level: intermediate

275: .seealso:  TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy()
276: J*/
277: typedef const char* TSTrajectoryType;
278: #define TSTRAJECTORYBASIC         "basic"
279: #define TSTRAJECTORYSINGLEFILE    "singlefile"
280: #define TSTRAJECTORYMEMORY        "memory"
281: #define TSTRAJECTORYVISUALIZATION "visualization"

283: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
284: PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
285: PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;

287: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
288: PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);

290: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
291: PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
292: PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
293: PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory,PetscViewer);
294: PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,TSTrajectoryType);
295: PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
296: PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
297: PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory,TS,PetscInt,PetscReal*,Vec,Vec);
298: PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory,TS,PetscReal,Vec*,Vec*);
299: PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory,PetscInt*);
300: PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory,Vec*,Vec*);
301: PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
302: PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[],PetscErrorCode (*)(TSTrajectory,TS));
303: PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
304: PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
305: PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory,PetscBool);
306: PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory,const char * const*);
307: PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
308: PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory,PetscBool);
309: PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory,PetscBool*);
310: PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory,PetscBool);
311: PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory,const char[]);
312: PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory,const char[]);
313: PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS,TSTrajectory*);

315: PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
316: PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
317: PETSC_EXTERN PetscErrorCode TSSetCostIntegrand(TS,PetscInt,Vec,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscBool,void*);
318: PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
319: PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec);

321: PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems*,TS);
322: PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
323: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
324: PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
325: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));

327: PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
328: PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
329: PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);

331: PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
332: PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
333: PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
334: PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
335: PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
336: PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);

338: PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS,PetscInt,Mat);
339: PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS,PetscInt*,Mat*);
340: PETSC_EXTERN PetscErrorCode TSForwardSetIntegralGradients(TS,PetscInt,Vec *);
341: PETSC_EXTERN PetscErrorCode TSForwardGetIntegralGradients(TS,PetscInt*,Vec **);
342: PETSC_EXTERN PetscErrorCode TSForwardSetRHSJacobianP(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,void*),void*);
343: PETSC_EXTERN PetscErrorCode TSForwardComputeRHSJacobianP(TS,PetscReal,Vec,Vec*);
344: PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
345: PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
346: PETSC_EXTERN PetscErrorCode TSForwardStep(TS);

348: PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS,PetscInt);
349: PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS,PetscInt*);
350: PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS,PetscReal);
351: PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS,PetscReal*);
352: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
353: PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS,TSExactFinalTimeOption*);

355: PETSC_EXTERN PETSC_DEPRECATED("Use TSSetTime[Step]")      PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
356: PETSC_EXTERN PETSC_DEPRECATED("Use TSSetMax{Steps|Time}") PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
357: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetMax{Steps|Time}") PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
358: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber")      PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
359: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber")      PetscErrorCode TSGetTotalSteps(TS,PetscInt*);

361: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
362: PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);

364: typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
365: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
366: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
367: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
368: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
369: PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
370: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS,PetscInt,PetscReal,Vec,void*);

372: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat *);
373: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);

375: PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
376: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
377: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);

379: PETSC_EXTERN PetscErrorCode TSStep(TS);
380: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*);
381: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
382: PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
383: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
384: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
385: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
386: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
387: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
388: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
389: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
390: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
391: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
392: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
393: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
394: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
395: PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
396: PETSC_EXTERN PetscErrorCode TSRollBack(TS);

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

400: PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
401: PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
402: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
403: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
404: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
405: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS,PetscInt*);
406: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS,PetscInt);

408: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
409: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
410: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
411: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
412: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
413: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
414: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);

416: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
417: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
418: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*);
419: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,TSForcingFunction,void*);

421: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
422: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
423: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
424: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
425: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
426: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);

428: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS,PetscReal,Vec,Vec,Vec,Vec,void*);
429: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat,void*);
430: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS,Vec,TSI2Function,void*);
431: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS,Vec*,TSI2Function*,void**);
432: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS,Mat,Mat,TSI2Jacobian,void*);
433: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS,Mat*,Mat*,TSI2Jacobian*,void**);

435: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS,const char[],IS);
436: PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS,const char[],IS*);
437: PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS,const char[],Vec,TSRHSFunction,void*);
438: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS,const char[],TS*);
439: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS,PetscInt*,TS*[]);

441: PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
442: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
443: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
444: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
445: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
446: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
447: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);

449: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
450: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
451: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
452: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode(*)(TS));
453: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
454: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
455: PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
456: PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
457: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
458: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
459: PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
460: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
461: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
462: PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
463: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
464: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
465: PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
466: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
467: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS,Vec,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
468: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
469: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
470: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
471: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);

473: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
474: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
475: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
476: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
477: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
478: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
479: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
480: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
481: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);

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

485: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
486: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
487: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
488: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
489: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS,PetscReal,Vec,Vec,Vec,Vec);
490: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat);
491: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);

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

495: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
496: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
497: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
498: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
499: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
500: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
501: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
502: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
503: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
504: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM,TSI2Function,void*);
505: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM,TSI2Function*,void**);
506: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM,TSI2Jacobian,void*);
507: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM,TSI2Jacobian*,void**);

509: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
510: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
511: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,TSForcingFunction,void*);
512: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,TSForcingFunction*,void**);
513: PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
514: PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
515: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **);

517: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
518: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
519: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);

521: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
522: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));

524: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
525: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
526: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
527: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);

529: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
530: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
531: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
532: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);

534: PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
535: PETSC_EXTERN PetscErrorCode DMPlexTSGetGradientDM(DM,PetscFV,DM*);

537: typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
538: typedef struct {
539:   Vec            ray;
540:   VecScatter     scatter;
541:   PetscViewer    viewer;
542:   TSMonitorLGCtx lgctx;
543: } TSMonitorDMDARayCtx;
544: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
545: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
546: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);

548: /* Dynamic creation and loading functions */
549: PETSC_EXTERN PetscFunctionList TSList;
550: PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
551: PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
552: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));

554: PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
555: PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
556: PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);

558: PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
559: PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
560: PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
561: PETSC_STATIC_INLINE PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}

563: #define TS_FILE_CLASSID 1211225

565: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
566: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);

568: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
569: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
570: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
571: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
572: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
573: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
574: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
575: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
576: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
577: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
578: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
579: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
580: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
581: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
582: PETSC_EXTERN PetscErrorCode TSMonitorError(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat *);

584: typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
585: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
586: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
587: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
588: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);

590: typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
591: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
592: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
593: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);

595: typedef struct _n_TSMonitorSPCtx* TSMonitorSPCtx;
596: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPCtx*);
597: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx*);
598: PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS,PetscInt,PetscReal,Vec,void*);

600: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
601: PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS,PetscReal);
602: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]);
603: /*J
604:    TSSSPType - string with the name of TSSSP scheme.

606:    Level: beginner

608: .seealso: TSSSPSetType(), TS
609: J*/
610: typedef const char* TSSSPType;
611: #define TSSSPRKS2  "rks2"
612: #define TSSSPRKS3  "rks3"
613: #define TSSSPRK104 "rk104"

615: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
616: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
617: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
618: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
619: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
620: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
621: PETSC_EXTERN PetscFunctionList TSSSPList;

623: /*S
624:    TSAdapt - Abstract object that manages time-step adaptivity

626:    Level: beginner

628: .seealso: TS, TSAdaptCreate(), TSAdaptType
629: S*/
630: typedef struct _p_TSAdapt *TSAdapt;

632: /*E
633:     TSAdaptType - String with the name of TSAdapt scheme.

635:    Level: beginner

637: .seealso: TSAdaptSetType(), TS
638: E*/
639: typedef const char *TSAdaptType;
640: #define TSADAPTNONE    "none"
641: #define TSADAPTBASIC   "basic"
642: #define TSADAPTDSP     "dsp"
643: #define TSADAPTCFL     "cfl"
644: #define TSADAPTGLEE    "glee"
645: #define TSADAPTHISTORY "history"

647: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
648: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
649: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
650: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
651: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
652: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
653: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt,TSAdaptType*);
654: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
655: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
656: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
657: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
658: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
659: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
660: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
661: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
662: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
663: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
664: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
665: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
666: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt,PetscBool);
667: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt,PetscReal,PetscReal);
668: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt,PetscReal*,PetscReal*);
669: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt,PetscReal,PetscReal);
670: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt,PetscReal*,PetscReal*);
671: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
672: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt,PetscReal*,PetscReal*);
673: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
674: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt,PetscInt n,PetscReal hist[],PetscBool);
675: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt,TSTrajectory,PetscBool);
676: PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt,PetscInt,PetscReal*,PetscReal*);
677: PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt,PetscInt);
678: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt,const char *);
679: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt,PetscReal,PetscReal,PetscReal);

681: /*S
682:    TSGLLEAdapt - Abstract object that manages time-step adaptivity

684:    Level: beginner

686:    Developer Notes:
687:    This functionality should be replaced by the TSAdapt.

689: .seealso: TSGLLE, TSGLLEAdaptCreate(), TSGLLEAdaptType
690: S*/
691: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;

693: /*J
694:     TSGLLEAdaptType - String with the name of TSGLLEAdapt scheme

696:    Level: beginner

698: .seealso: TSGLLEAdaptSetType(), TS
699: J*/
700: typedef const char *TSGLLEAdaptType;
701: #define TSGLLEADAPT_NONE "none"
702: #define TSGLLEADAPT_SIZE "size"
703: #define TSGLLEADAPT_BOTH "both"

705: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[],PetscErrorCode (*)(TSGLLEAdapt));
706: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
707: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
708: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm,TSGLLEAdapt*);
709: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt,TSGLLEAdaptType);
710: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt,const char[]);
711: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
712: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt,PetscViewer);
713: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems*,TSGLLEAdapt);
714: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt*);

716: /*J
717:     TSGLLEAcceptType - String with the name of TSGLLEAccept scheme

719:    Level: beginner

721: .seealso: TSGLLESetAcceptType(), TS
722: J*/
723: typedef const char *TSGLLEAcceptType;
724: #define TSGLLEACCEPT_ALWAYS "always"

726: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
727: PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[],TSGLLEAcceptFunction);

729: /*J
730:   TSGLLEType - family of time integration method within the General Linear class

732:   Level: beginner

734: .seealso: TSGLLESetType(), TSGLLERegister()
735: J*/
736: typedef const char* TSGLLEType;
737: #define TSGLLE_IRKS   "irks"

739: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[],PetscErrorCode(*)(TS));
740: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
741: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
742: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS,TSGLLEType);
743: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS,TSGLLEAdapt*);
744: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS,TSGLLEAcceptType);

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

749:    Level: beginner

751: .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
752: J*/
753: #define TSEIMEXType   char*

755: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
756: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
757: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);

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

762:    Level: beginner

764: .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
765: J*/
766: typedef const char* TSRKType;
767: #define TSRK1FE   "1fe"
768: #define TSRK2A    "2a"
769: #define TSRK3     "3"
770: #define TSRK3BS   "3bs"
771: #define TSRK4     "4"
772: #define TSRK5F    "5f"
773: #define TSRK5DP   "5dp"
774: #define TSRK5BS   "5bs"

776: PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
777: PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
778: PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
779: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
780: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
781: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
782: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);

784: /*J
785:     TSGLEEType - String with the name of a General Linear with Error Estimation method.

787:    Level: beginner

789: .seealso: TSGLEESetType(), TS, TSGLEE, TSGLEERegister()
790: J*/
791: typedef const char* TSGLEEType;
792: #define TSGLEEi1      "BE1"
793: #define TSGLEE23      "23"
794: #define TSGLEE24      "24"
795: #define TSGLEE25I     "25i"
796: #define TSGLEE35      "35"
797: #define TSGLEEEXRK2A  "exrk2a"
798: #define TSGLEERK32G1  "rk32g1"
799: #define TSGLEERK285EX "rk285ex"
800: /*J
801:     TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method.

803:    Level: beginner

805: .seealso: TSGLEESetMode(), TS, TSGLEE, TSGLEERegister()
806: J*/
807: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType*);
808: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts,TSGLEEType);
809: PETSC_EXTERN PetscErrorCode TSGLEERegister(TSGLEEType,PetscInt,PetscInt,PetscInt,PetscReal,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
810: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
811: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
812: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);

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

817:    Level: beginner

819: .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
820: J*/
821: typedef const char* TSARKIMEXType;
822: #define TSARKIMEX1BEE   "1bee"
823: #define TSARKIMEXA2     "a2"
824: #define TSARKIMEXL2     "l2"
825: #define TSARKIMEXARS122 "ars122"
826: #define TSARKIMEX2C     "2c"
827: #define TSARKIMEX2D     "2d"
828: #define TSARKIMEX2E     "2e"
829: #define TSARKIMEXPRSSP2 "prssp2"
830: #define TSARKIMEX3      "3"
831: #define TSARKIMEXBPR3   "bpr3"
832: #define TSARKIMEXARS443 "ars443"
833: #define TSARKIMEX4      "4"
834: #define TSARKIMEX5      "5"
835: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
836: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
837: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
838: 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[]);
839: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
840: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
841: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);

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

846:    Level: beginner

848: .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
849: J*/
850: typedef const char* TSRosWType;
851: #define TSROSW2M          "2m"
852: #define TSROSW2P          "2p"
853: #define TSROSWRA3PW       "ra3pw"
854: #define TSROSWRA34PW2     "ra34pw2"
855: #define TSROSWRODAS3      "rodas3"
856: #define TSROSWSANDU3      "sandu3"
857: #define TSROSWASSP3P3S1C  "assp3p3s1c"
858: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
859: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
860: #define TSROSWARK3        "ark3"
861: #define TSROSWTHETA1      "theta1"
862: #define TSROSWTHETA2      "theta2"
863: #define TSROSWGRK4T       "grk4t"
864: #define TSROSWSHAMP4      "shamp4"
865: #define TSROSWVELDD4      "veldd4"
866: #define TSROSW4L          "4l"

868: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS,TSRosWType*);
869: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS,TSRosWType);
870: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
871: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
872: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
873: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
874: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
875: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);

877: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS,PetscInt);
878: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS,PetscInt*);

880: /*J
881:   TSBasicSymplecticType - String with the name of a basic symplectic integration method.

883:   Level: beginner

885:   .seealso: TSBasicSymplecticSetType(), TS, TSBASICSYMPLECTIC, TSBasicSymplecticRegister()
886: J*/
887: typedef const char* TSBasicSymplecticType;
888: #define TSBASICSYMPLECTICSIEULER   "1"
889: #define TSBASICSYMPLECTICVELVERLET "2"
890: #define TSBASICSYMPLECTIC3         "3"
891: #define TSBASICSYMPLECTIC4         "4"
892: PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS,TSBasicSymplecticType);
893: PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS,TSBasicSymplecticType*);
894: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType,PetscInt,PetscInt,PetscReal[],PetscReal[]);
895: PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
896: PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
897: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);

899: /*
900:        PETSc interface to Sundials
901: */
902: #ifdef PETSC_HAVE_SUNDIALS
903: typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
904: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
905: typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
906: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
907: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
908: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
909: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
910: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
911: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
912: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
913: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
914: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
915: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
916: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
917: PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
918: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
919: #endif

921: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
922: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
923: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
924: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);

926: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
927: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
928: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);

930: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS,PetscReal);
931: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS,PetscReal,PetscReal,PetscReal,PetscReal);
932: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*);

934: PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
935: PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);

937: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
938: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);

940: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS,PetscBool*);
941: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS,PetscBool*);
942: #endif