Actual source code: petscts.h

petsc-3.9.4 2018-09-11
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 TSALPHA2          "alpha2"
 38: #define TSGLLE            "glle"
 39: #define TSGLEE            "glee"
 40: #define TSSSP             "ssp"
 41: #define TSARKIMEX         "arkimex"
 42: #define TSROSW            "rosw"
 43: #define TSEIMEX           "eimex"
 44: #define TSMIMEX           "mimex"
 45: #define TSBDF             "bdf"
 46: #define TSRADAU5          "radau5"

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

 51:    Level: beginner

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

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

 60:    Level: beginner

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

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

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

 88: /*E
 89:    TSConvergedReason - reason a TS method has converged or not

 91:    Level: beginner

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

 95:    Each reason has its own manual page.

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

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

115:    Level: beginner

117: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
118: M*/

120: /*MC
121:    TS_CONVERGED_TIME - the final time was reached

123:    Level: beginner

125: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxTime(), TSGetMaxTime(), TSGetSolveTime()
126: M*/

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

131:    Level: beginner

133: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxSteps(), TSGetMaxSteps()
134: M*/

136: /*MC
137:    TS_CONVERGED_USER - user requested termination

139:    Level: beginner

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

144: /*MC
145:    TS_CONVERGED_EVENT - user requested termination on event detection

147:    Level: beginner

149: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
150: M*/

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

155:    Level: beginner

157:    Options Database:
158: .   -ts_pseudo_frtol <rtol>

160: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FATOL
161: M*/

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

166:    Level: beginner

168:    Options Database:
169: .   -ts_pseudo_fatol <atol>

171: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FRTOL
172: M*/

174: /*MC
175:    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed

177:    Level: beginner

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

181: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
182: M*/

184: /*MC
185:    TS_DIVERGED_STEP_REJECTED - too many steps were rejected

187:    Level: beginner

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

191: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
192: M*/

194: /*E
195:    TSExactFinalTimeOption - option for handling of final time step

197:    Level: beginner

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

201: $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
202: $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
203: $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time

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

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

211: /* Logging support */
212: PETSC_EXTERN PetscClassId TS_CLASSID;
213: PETSC_EXTERN PetscClassId DMTS_CLASSID;
214: PETSC_EXTERN PetscClassId TSADAPT_CLASSID;

216: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
217: PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);

219: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
220: PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
221: PETSC_EXTERN PetscErrorCode TSDestroy(TS*);

223: PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
224: PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
225: PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
226: PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
227: PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
228: PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);

230: PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
231: PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
232: PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
233: PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
234: PETSC_EXTERN PetscErrorCode TSSetUp(TS);
235: PETSC_EXTERN PetscErrorCode TSReset(TS);

237: PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
238: PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);

240: PETSC_EXTERN PetscErrorCode TS2SetSolution(TS,Vec,Vec);
241: PETSC_EXTERN PetscErrorCode TS2GetSolution(TS,Vec*,Vec*);

243: PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS,PetscInt*,Vec*);
244: PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS,Vec*);
245: PETSC_EXTERN PetscErrorCode TSGetTimeError(TS,PetscInt,Vec*);
246: PETSC_EXTERN PetscErrorCode TSSetTimeError(TS,Vec);

248: PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
249: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS,PetscReal,Vec,Mat);
250: PETSC_EXTERN PetscErrorCode TSComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
251: PETSC_EXTERN PetscErrorCode TSComputeDRDYFunction(TS,PetscReal,Vec,Vec*);

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

256:    Level: advanced

258:   Concepts: ODE solvers, trajectory

260: .seealso:  TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectorySetType(), TSTrajectoryDestroy()
261: S*/
262: typedef struct _p_TSTrajectory* TSTrajectory;

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

267:    Level: intermediate

269: .seealso:  TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy()
270: J*/
271: typedef const char* TSTrajectoryType;
272: #define TSTRAJECTORYBASIC         "basic"
273: #define TSTRAJECTORYSINGLEFILE    "singlefile"
274: #define TSTRAJECTORYMEMORY        "memory"
275: #define TSTRAJECTORYVISUALIZATION "visualization"

277: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
278: PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
279: PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;

281: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);

283: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
284: PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
285: PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory,PetscViewer);
286: PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,const TSTrajectoryType);
287: PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
288: PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
289: PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
290: PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[],PetscErrorCode (*)(TSTrajectory,TS));
291: PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
292: PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
293: PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory,PetscBool);
294: PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory,const char * const*);
295: PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
296: PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory,PetscBool);
297: PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory,const char[]);
298: PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory,const char[]);
299: PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS,TSTrajectory*);

301: PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
302: PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
303: 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*);
304: PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
305: PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec);

307: PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems*,TS);
308: PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
309: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
310: PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
311: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));

313: PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
314: PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
315: PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);

317: PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
318: PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
319: PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
320: PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
321: PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
322: PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);

324: PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS,PetscInt,Mat);
325: PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS,PetscInt*,Mat*);
326: PETSC_EXTERN PetscErrorCode TSForwardSetIntegralGradients(TS,PetscInt,Vec *);
327: PETSC_EXTERN PetscErrorCode TSForwardGetIntegralGradients(TS,PetscInt*,Vec **);
328: PETSC_EXTERN PetscErrorCode TSForwardSetRHSJacobianP(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,void*),void*);
329: PETSC_EXTERN PetscErrorCode TSForwardComputeRHSJacobianP(TS,PetscReal,Vec,Vec*);
330: PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
331: PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
332: PETSC_EXTERN PetscErrorCode TSForwardStep(TS);

334: PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS,PetscInt);
335: PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS,PetscInt*);
336: PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS,PetscReal);
337: PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS,PetscReal*);
338: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
339: PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS,TSExactFinalTimeOption*);

341: PETSC_EXTERN PETSC_DEPRECATED("Use TSSetTime[Step]")      PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
342: PETSC_EXTERN PETSC_DEPRECATED("Use TSSetMax{Steps|Time}") PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
343: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetMax{Steps|Time}") PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
344: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber")      PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
345: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber")      PetscErrorCode TSGetTotalSteps(TS,PetscInt*);

347: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);

349: typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
350: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
351: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
352: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
353: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
354: PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
355: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS,PetscInt,PetscReal,Vec,void*);

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

360: PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
361: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
362: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);

364: PETSC_EXTERN PetscErrorCode TSStep(TS);
365: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*);
366: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
367: PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
368: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
369: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
370: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
371: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
372: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
373: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
374: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
375: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
376: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
377: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
378: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
379: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
380: PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
381: PETSC_EXTERN PetscErrorCode TSRollBack(TS);

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

385: PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
386: PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
387: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
388: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
389: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
390: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS,PetscInt*);
391: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS,PetscInt);

393: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
394: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
395: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
396: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
397: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
398: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
399: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);

401: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
402: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
403: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*);
404: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,TSForcingFunction,void*);

406: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
407: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
408: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
409: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
410: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
411: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);

413: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS,PetscReal,Vec,Vec,Vec,Vec,void*);
414: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat,void*);
415: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS,Vec,TSI2Function,void*);
416: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS,Vec*,TSI2Function*,void**);
417: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS,Mat,Mat,TSI2Jacobian,void*);
418: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS,Mat*,Mat*,TSI2Jacobian*,void**);

420: PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
421: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
422: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
423: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
424: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
425: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
426: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);

428: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
429: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
430: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
431: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode(*)(TS));
432: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
433: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
434: PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
435: PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
436: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
437: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
438: PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
439: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
440: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
441: PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
442: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
443: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
444: PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
445: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
446: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS,Vec,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
447: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
448: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
449: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
450: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);

452: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
453: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
454: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
455: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
456: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
457: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
458: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
459: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
460: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);

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

464: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
465: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
466: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
467: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
468: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS,PetscReal,Vec,Vec,Vec,Vec);
469: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat);
470: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);

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

474: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
475: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
476: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
477: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
478: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
479: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
480: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
481: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
482: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
483: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM,TSI2Function,void*);
484: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM,TSI2Function*,void**);
485: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM,TSI2Jacobian,void*);
486: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM,TSI2Jacobian*,void**);

488: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
489: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
490: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,TSForcingFunction,void*);
491: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,TSForcingFunction*,void**);
492: PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
493: PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
494: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **);

496: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
497: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
498: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);

500: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
501: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));

503: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
504: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
505: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
506: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);

508: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
509: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
510: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
511: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);

513: PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
514: PETSC_EXTERN PetscErrorCode DMPlexTSGetGradientDM(DM,PetscFV,DM*);

516: typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
517: typedef struct {
518:   Vec            ray;
519:   VecScatter     scatter;
520:   PetscViewer    viewer;
521:   TSMonitorLGCtx lgctx;
522: } TSMonitorDMDARayCtx;
523: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
524: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
525: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);

527: /* Dynamic creation and loading functions */
528: PETSC_EXTERN PetscFunctionList TSList;
529: PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
530: PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
531: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));

533: PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
534: PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
535: PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);

537: PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
538: PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
539: PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
540: PETSC_STATIC_INLINE PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}

542: #define TS_FILE_CLASSID 1211225

544: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
545: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);

547: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
548: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
549: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
550: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
551: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
552: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
553: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
554: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
555: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
556: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
557: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
558: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
559: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
560: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
561: PETSC_EXTERN PetscErrorCode TSMonitorError(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat *);

563: typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
564: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
565: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
566: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
567: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);

569: typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
570: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
571: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
572: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);

574: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
575: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]);
576: /*J
577:    TSSSPType - string with the name of TSSSP scheme.

579:    Level: beginner

581: .seealso: TSSSPSetType(), TS
582: J*/
583: typedef const char* TSSSPType;
584: #define TSSSPRKS2  "rks2"
585: #define TSSSPRKS3  "rks3"
586: #define TSSSPRK104 "rk104"

588: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
589: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
590: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
591: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
592: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
593: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
594: PETSC_EXTERN PetscFunctionList TSSSPList;

596: /*S
597:    TSAdapt - Abstract object that manages time-step adaptivity

599:    Level: beginner

601: .seealso: TS, TSAdaptCreate(), TSAdaptType
602: S*/
603: typedef struct _p_TSAdapt *TSAdapt;

605: /*E
606:     TSAdaptType - String with the name of TSAdapt scheme.

608:    Level: beginner

610: .seealso: TSAdaptSetType(), TS
611: E*/
612: typedef const char *TSAdaptType;
613: #define TSADAPTNONE  "none"
614: #define TSADAPTBASIC "basic"
615: #define TSADAPTDSP   "dsp"
616: #define TSADAPTCFL   "cfl"
617: #define TSADAPTGLEE  "glee"

619: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
620: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
621: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
622: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
623: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
624: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
625: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt,TSAdaptType*);
626: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
627: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
628: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
629: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
630: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
631: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
632: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
633: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
634: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
635: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
636: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
637: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
638: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt,PetscBool);
639: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt,PetscReal,PetscReal);
640: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt,PetscReal*,PetscReal*);
641: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt,PetscReal,PetscReal);
642: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt,PetscReal*,PetscReal*);
643: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
644: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt,PetscReal*,PetscReal*);
645: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
646: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt,const char *);
647: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt,PetscReal,PetscReal,PetscReal);

649: /*S
650:    TSGLLEAdapt - Abstract object that manages time-step adaptivity

652:    Level: beginner

654:    Developer Notes:
655:    This functionality should be replaced by the TSAdapt.

657: .seealso: TSGLLE, TSGLLEAdaptCreate(), TSGLLEAdaptType
658: S*/
659: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;

661: /*J
662:     TSGLLEAdaptType - String with the name of TSGLLEAdapt scheme

664:    Level: beginner

666: .seealso: TSGLLEAdaptSetType(), TS
667: J*/
668: typedef const char *TSGLLEAdaptType;
669: #define TSGLLEADAPT_NONE "none"
670: #define TSGLLEADAPT_SIZE "size"
671: #define TSGLLEADAPT_BOTH "both"

673: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[],PetscErrorCode (*)(TSGLLEAdapt));
674: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
675: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
676: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm,TSGLLEAdapt*);
677: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt,TSGLLEAdaptType);
678: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt,const char[]);
679: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
680: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt,PetscViewer);
681: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems*,TSGLLEAdapt);
682: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt*);

684: /*J
685:     TSGLLEAcceptType - String with the name of TSGLLEAccept scheme

687:    Level: beginner

689: .seealso: TSGLLESetAcceptType(), TS
690: J*/
691: typedef const char *TSGLLEAcceptType;
692: #define TSGLLEACCEPT_ALWAYS "always"

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

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

700:   Level: beginner

702: .seealso: TSGLLESetType(), TSGLLERegister()
703: J*/
704: typedef const char* TSGLLEType;
705: #define TSGLLE_IRKS   "irks"

707: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[],PetscErrorCode(*)(TS));
708: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
709: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
710: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS,TSGLLEType);
711: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS,TSGLLEAdapt*);
712: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS,TSGLLEAcceptType);

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

717:    Level: beginner

719: .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
720: J*/
721: #define TSEIMEXType   char*

723: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
724: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
725: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);

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

730:    Level: beginner

732: .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
733: J*/
734: typedef const char* TSRKType;
735: #define TSRK1FE   "1fe"
736: #define TSRK2A    "2a"
737: #define TSRK3     "3"
738: #define TSRK3BS   "3bs"
739: #define TSRK4     "4"
740: #define TSRK5F    "5f"
741: #define TSRK5DP   "5dp"
742: #define TSRK5BS   "5bs"

744: PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
745: PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
746: PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
747: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
748: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
749: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
750: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);

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

755:    Level: beginner

757: .seealso: TSGLEESetType(), TS, TSGLEE, TSGLEERegister()
758: J*/
759: typedef const char* TSGLEEType;
760: #define TSGLEEi1      "BE1"
761: #define TSGLEE23      "23"
762: #define TSGLEE24      "24"
763: #define TSGLEE25I     "25i"
764: #define TSGLEE35      "35"
765: #define TSGLEEEXRK2A  "exrk2a"
766: #define TSGLEERK32G1  "rk32g1"
767: #define TSGLEERK285EX "rk285ex"
768: /*J
769:     TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method.

771:    Level: beginner

773: .seealso: TSGLEESetMode(), TS, TSGLEE, TSGLEERegister()
774: J*/
775: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType*);
776: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts,TSGLEEType);
777: 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[]);
778: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
779: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
780: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);

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

785:    Level: beginner

787: .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
788: J*/
789: typedef const char* TSARKIMEXType;
790: #define TSARKIMEX1BEE   "1bee"
791: #define TSARKIMEXA2     "a2"
792: #define TSARKIMEXL2     "l2"
793: #define TSARKIMEXARS122 "ars122"
794: #define TSARKIMEX2C     "2c"
795: #define TSARKIMEX2D     "2d"
796: #define TSARKIMEX2E     "2e"
797: #define TSARKIMEXPRSSP2 "prssp2"
798: #define TSARKIMEX3      "3"
799: #define TSARKIMEXBPR3   "bpr3"
800: #define TSARKIMEXARS443 "ars443"
801: #define TSARKIMEX4      "4"
802: #define TSARKIMEX5      "5"
803: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
804: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
805: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
806: 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[]);
807: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
808: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
809: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);

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

814:    Level: beginner

816: .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
817: J*/
818: typedef const char* TSRosWType;
819: #define TSROSW2M          "2m"
820: #define TSROSW2P          "2p"
821: #define TSROSWRA3PW       "ra3pw"
822: #define TSROSWRA34PW2     "ra34pw2"
823: #define TSROSWRODAS3      "rodas3"
824: #define TSROSWSANDU3      "sandu3"
825: #define TSROSWASSP3P3S1C  "assp3p3s1c"
826: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
827: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
828: #define TSROSWARK3        "ark3"
829: #define TSROSWTHETA1      "theta1"
830: #define TSROSWTHETA2      "theta2"
831: #define TSROSWGRK4T       "grk4t"
832: #define TSROSWSHAMP4      "shamp4"
833: #define TSROSWVELDD4      "veldd4"
834: #define TSROSW4L          "4l"

836: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
837: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
838: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
839: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
840: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
841: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
842: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
843: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);

845: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS,PetscInt);
846: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS,PetscInt*);

848: /*
849:        PETSc interface to Sundials
850: */
851: #ifdef PETSC_HAVE_SUNDIALS
852: typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
853: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
854: typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
855: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
856: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
857: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
858: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
859: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
860: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
861: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
862: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
863: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
864: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
865: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
866: PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
867: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
868: #endif

870: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
871: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
872: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
873: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);

875: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
876: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
877: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);

879: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS,PetscReal);
880: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS,PetscReal,PetscReal,PetscReal,PetscReal);
881: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*);

883: PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
884: PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);

886: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
887: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);

889: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS,PetscBool*);
890: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS,PetscBool*);
891: #endif