Actual source code: petscts.h

petsc-master 2019-12-13
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>
  8:  #include <petscconvest.h>

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

 13:    Level: beginner

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

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

 22:    Level: beginner

 24: .seealso: TSSetType(), TS, TSRegister()
 25: J*/
 26: typedef const char* TSType;
 27: #define TSEULER           "euler"
 28: #define TSBEULER          "beuler"
 29: #define TSBASICSYMPLECTIC "basicsymplectic"
 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"
 47: #define TSMPRK            "mprk"

 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:   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
113:   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
114: } TSConvergedReason;
115: PETSC_EXTERN const char *const*TSConvergedReasons;

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

120:    Level: beginner

122: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
123: M*/

125: /*MC
126:    TS_CONVERGED_TIME - the final time was reached

128:    Level: beginner

130: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxTime(), TSGetMaxTime(), TSGetSolveTime()
131: M*/

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

136:    Level: beginner

138: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxSteps(), TSGetMaxSteps()
139: M*/

141: /*MC
142:    TS_CONVERGED_USER - user requested termination

144:    Level: beginner

146: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
147: M*/

149: /*MC
150:    TS_CONVERGED_EVENT - user requested termination on event detection

152:    Level: beginner

154: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
155: M*/

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

160:    Level: beginner

162:    Options Database:
163: .   -ts_pseudo_frtol <rtol>

165: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FATOL
166: M*/

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

171:    Level: beginner

173:    Options Database:
174: .   -ts_pseudo_fatol <atol>

176: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FRTOL
177: M*/

179: /*MC
180:    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed

182:    Level: beginner

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

187: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
188: M*/

190: /*MC
191:    TS_DIVERGED_STEP_REJECTED - too many steps were rejected

193:    Level: beginner

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

198: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
199: M*/

201: /*E
202:    TSExactFinalTimeOption - option for handling of final time step

204:    Level: beginner

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

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

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

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

219: /* Logging support */
220: PETSC_EXTERN PetscClassId TS_CLASSID;
221: PETSC_EXTERN PetscClassId DMTS_CLASSID;
222: PETSC_EXTERN PetscClassId TSADAPT_CLASSID;

224: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
225: PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);

227: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
228: PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
229: PETSC_EXTERN PetscErrorCode TSDestroy(TS*);

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

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

245: PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
246: PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);

248: PETSC_EXTERN PetscErrorCode TS2SetSolution(TS,Vec,Vec);
249: PETSC_EXTERN PetscErrorCode TS2GetSolution(TS,Vec*,Vec*);

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

256: PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
257: PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS,Mat*,PetscErrorCode(**)(TS,PetscReal,Vec,Mat,void*),void**);
258: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS,PetscReal,Vec,Mat);
259: PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void*);
260: PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS,PetscReal,Vec,Vec,PetscReal,Mat,PetscBool);
261: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
262: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSComputeDRDUFunction(TS,PetscReal,Vec,Vec*);
263: PETSC_EXTERN PetscErrorCode TSSetIHessianProduct(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
264:                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
265:                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
266:                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
267:                                                     void*);
268: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
269: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
270: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
271: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
272: PETSC_EXTERN PetscErrorCode TSSetRHSHessianProduct(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
273:                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
274:                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
275:                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
276:                                                     void*);
277: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
278: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
279: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
280: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
281: PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS,PetscInt,Vec*,Vec*,Vec);
282: PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS,PetscInt*,Vec**,Vec**,Vec*);
283: PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS,Vec,Mat,Mat);

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

288:    Level: advanced

290: .seealso:  TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectorySetType(), TSTrajectoryDestroy(), TSTrajectoryReset()
291: S*/
292: typedef struct _p_TSTrajectory* TSTrajectory;

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

297:    Level: intermediate

299: .seealso:  TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy()
300: J*/
301: typedef const char* TSTrajectoryType;
302: #define TSTRAJECTORYBASIC         "basic"
303: #define TSTRAJECTORYSINGLEFILE    "singlefile"
304: #define TSTRAJECTORYMEMORY        "memory"
305: #define TSTRAJECTORYVISUALIZATION "visualization"

307: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
308: PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
309: PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;

311: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
312: PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);

314: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
315: PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
316: PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
317: PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory,PetscViewer);
318: PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,TSTrajectoryType);
319: PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory,TS,TSTrajectoryType*);
320: PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
321: PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
322: PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory,TS,PetscInt,PetscReal*,Vec,Vec);
323: PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory,TS,PetscReal,Vec*,Vec*);
324: PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory,PetscInt*);
325: PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory,Vec*,Vec*);
326: PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
327: PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[],PetscErrorCode (*)(TSTrajectory,TS));
328: PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
329: PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
330: PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory,PetscBool);
331: PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory,PetscBool);
332: PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory,const char * const*);
333: PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
334: PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory,PetscBool);
335: PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory,PetscBool*);
336: PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory,PetscBool);
337: PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory,const char[]);
338: PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory,const char[]);
339: PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS,TSTrajectory*);

341: PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
342: PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
343: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSCreateQuadratureTS() then set up the sub-TS (since version 3.12)") 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*);
344: PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
345: PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec);
346: PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS,PetscBool,TS*);
347: PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS,PetscBool*,TS*);

349: PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems*,TS);
350: PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
351: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
352: PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
353: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));

355: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetRHSJacobianP()")     PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
356: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSComputeRHSJacobianP()") PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
357: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
358: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
359: PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
360: PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);

362: PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
363: PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
364: PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
365: PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
366: PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS,Mat);
367: PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);

369: PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS,PetscInt,Mat);
370: PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS,PetscInt*,Mat*);
371: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSCreateQuadratureTS() and TSForwardSetSensitivities() (since version 3.12)") PetscErrorCode TSForwardSetIntegralGradients(TS,PetscInt,Vec *);
372: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSForwardGetSensitivities()")                          PetscErrorCode TSForwardGetIntegralGradients(TS,PetscInt*,Vec **);
373: PETSC_EXTERN PetscErrorCode TSForwardSetRHSJacobianP(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,void*),void*);
374: PETSC_EXTERN PetscErrorCode TSForwardComputeRHSJacobianP(TS,PetscReal,Vec,Vec*);
375: PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
376: PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
377: PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
378: PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
379: PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS,Mat);
380: PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS,PetscInt*,Mat**);

382: PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS,PetscInt);
383: PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS,PetscInt*);
384: PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS,PetscReal);
385: PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS,PetscReal*);
386: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
387: PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS,TSExactFinalTimeOption*);

389: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetTime[Step]() (since version 3.8)")      PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
390: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
391: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
392: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)")      PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
393: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)")      PetscErrorCode TSGetTotalSteps(TS,PetscInt*);

395: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
396: PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);

398: typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
399: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
400: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
401: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
402: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
403: PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
404: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS,PetscInt,PetscReal,Vec,void*);

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

409: PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
410: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
411: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);

413: PETSC_EXTERN PetscErrorCode TSStep(TS);
414: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*);
415: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
416: PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
417: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
418: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
419: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
420: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
421: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
422: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
423: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
424: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
425: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
426: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
427: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
428: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
429: PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
430: PETSC_EXTERN PetscErrorCode TSRollBack(TS);

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

434: PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
435: PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
436: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
437: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
438: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
439: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS,PetscInt*);
440: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS,PetscInt);

442: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
443: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
444: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS,PetscReal,Vec,Mat,void*);
445: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
446: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
447: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
448: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
449: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);

451: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
452: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
453: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*);
454: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,TSForcingFunction,void*);

456: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
457: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
458: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
459: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
460: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
461: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);

463: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS,PetscReal,Vec,Vec,Vec,Vec,void*);
464: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat,void*);
465: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS,Vec,TSI2Function,void*);
466: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS,Vec*,TSI2Function*,void**);
467: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS,Mat,Mat,TSI2Jacobian,void*);
468: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS,Mat*,Mat*,TSI2Jacobian*,void**);

470: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS,const char[],IS);
471: PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS,const char[],IS*);
472: PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS,const char[],Vec,TSRHSFunction,void*);
473: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS,const char[],TS*);
474: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS,PetscInt*,TS*[]);
475: PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
476: PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool*);

478: PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
479: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
480: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
481: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
482: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
483: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
484: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);

486: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
487: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
488: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
489: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode(*)(TS));
490: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
491: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
492: PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
493: PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
494: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
495: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
496: PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
497: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
498: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
499: PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
500: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
501: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
502: PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
503: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
504: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS,Vec,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
505: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
506: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
507: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
508: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);

510: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
511: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
512: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
513: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
514: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
515: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
516: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
517: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
518: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);

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

522: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
523: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
524: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
525: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
526: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS,PetscReal,Vec,Vec,Vec,Vec);
527: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat);
528: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);

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

532: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
533: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
534: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
535: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
536: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
537: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
538: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
539: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
540: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
541: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM,TSI2Function,void*);
542: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM,TSI2Function*,void**);
543: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM,TSI2Jacobian,void*);
544: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM,TSI2Jacobian*,void**);

546: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
547: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
548: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,TSForcingFunction,void*);
549: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,TSForcingFunction*,void**);
550: PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
551: PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
552: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **);

554: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
555: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
556: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);

558: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
559: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));

561: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
562: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
563: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
564: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);

566: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
567: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
568: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
569: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);

571: PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
572: PETSC_EXTERN PetscErrorCode DMPlexTSGetGradientDM(DM,PetscFV,DM*);

574: typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
575: typedef struct {
576:   Vec            ray;
577:   VecScatter     scatter;
578:   PetscViewer    viewer;
579:   TSMonitorLGCtx lgctx;
580: } TSMonitorDMDARayCtx;
581: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
582: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
583: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);

585: /* Dynamic creation and loading functions */
586: PETSC_EXTERN PetscFunctionList TSList;
587: PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
588: PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
589: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));

591: PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
592: PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
593: PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);

595: PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
596: PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
597: PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS,PetscObject,const char[]);
598: PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory,PetscObject,const char[]);

600: #define TS_FILE_CLASSID 1211225

602: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
603: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);

605: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
606: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
607: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
608: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
609: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
610: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
611: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
612: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
613: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
614: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
615: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
616: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
617: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
618: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
619: PETSC_EXTERN PetscErrorCode TSMonitorError(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat *);

621: typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
622: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
623: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
624: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
625: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);

627: typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
628: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
629: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
630: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);

632: typedef struct _n_TSMonitorSPCtx* TSMonitorSPCtx;
633: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPCtx*);
634: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx*);
635: PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS,PetscInt,PetscReal,Vec,void*);

637: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
638: PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS,PetscReal);
639: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]);
640: /*J
641:    TSSSPType - string with the name of TSSSP scheme.

643:    Level: beginner

645: .seealso: TSSSPSetType(), TS
646: J*/
647: typedef const char* TSSSPType;
648: #define TSSSPRKS2  "rks2"
649: #define TSSSPRKS3  "rks3"
650: #define TSSSPRK104 "rk104"

652: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
653: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
654: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
655: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
656: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
657: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
658: PETSC_EXTERN PetscFunctionList TSSSPList;

660: /*S
661:    TSAdapt - Abstract object that manages time-step adaptivity

663:    Level: beginner

665: .seealso: TS, TSAdaptCreate(), TSAdaptType
666: S*/
667: typedef struct _p_TSAdapt *TSAdapt;

669: /*E
670:     TSAdaptType - String with the name of TSAdapt scheme.

672:    Level: beginner

674: .seealso: TSAdaptSetType(), TS
675: E*/
676: typedef const char *TSAdaptType;
677: #define TSADAPTNONE    "none"
678: #define TSADAPTBASIC   "basic"
679: #define TSADAPTDSP     "dsp"
680: #define TSADAPTCFL     "cfl"
681: #define TSADAPTGLEE    "glee"
682: #define TSADAPTHISTORY "history"

684: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
685: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
686: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
687: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
688: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
689: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
690: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt,TSAdaptType*);
691: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
692: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
693: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
694: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
695: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
696: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
697: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
698: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
699: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
700: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
701: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
702: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
703: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt,PetscBool);
704: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt,PetscReal,PetscReal);
705: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt,PetscReal*,PetscReal*);
706: PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt,PetscReal);
707: PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt,PetscReal*);
708: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt,PetscReal,PetscReal);
709: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt,PetscReal*,PetscReal*);
710: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
711: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt,PetscReal*,PetscReal*);
712: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
713: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt,PetscInt n,PetscReal hist[],PetscBool);
714: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt,TSTrajectory,PetscBool);
715: PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt,PetscInt,PetscReal*,PetscReal*);
716: PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt,PetscInt);
717: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt,const char *);
718: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt,PetscReal,PetscReal,PetscReal);

720: /*S
721:    TSGLLEAdapt - Abstract object that manages time-step adaptivity

723:    Level: beginner

725:    Developer Notes:
726:    This functionality should be replaced by the TSAdapt.

728: .seealso: TSGLLE, TSGLLEAdaptCreate(), TSGLLEAdaptType
729: S*/
730: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;

732: /*J
733:     TSGLLEAdaptType - String with the name of TSGLLEAdapt scheme

735:    Level: beginner

737: .seealso: TSGLLEAdaptSetType(), TS
738: J*/
739: typedef const char *TSGLLEAdaptType;
740: #define TSGLLEADAPT_NONE "none"
741: #define TSGLLEADAPT_SIZE "size"
742: #define TSGLLEADAPT_BOTH "both"

744: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[],PetscErrorCode (*)(TSGLLEAdapt));
745: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
746: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
747: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm,TSGLLEAdapt*);
748: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt,TSGLLEAdaptType);
749: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt,const char[]);
750: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
751: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt,PetscViewer);
752: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems*,TSGLLEAdapt);
753: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt*);

755: /*J
756:     TSGLLEAcceptType - String with the name of TSGLLEAccept scheme

758:    Level: beginner

760: .seealso: TSGLLESetAcceptType(), TS
761: J*/
762: typedef const char *TSGLLEAcceptType;
763: #define TSGLLEACCEPT_ALWAYS "always"

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

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

771:   Level: beginner

773: .seealso: TSGLLESetType(), TSGLLERegister()
774: J*/
775: typedef const char* TSGLLEType;
776: #define TSGLLE_IRKS   "irks"

778: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[],PetscErrorCode(*)(TS));
779: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
780: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
781: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS,TSGLLEType);
782: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS,TSGLLEAdapt*);
783: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS,TSGLLEAcceptType);

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

788:    Level: beginner

790: .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
791: J*/
792: #define TSEIMEXType   char*

794: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
795: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
796: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);

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

801:    Level: beginner

803: .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
804: J*/
805: typedef const char* TSRKType;
806: #define TSRK1FE   "1fe"
807: #define TSRK2A    "2a"
808: #define TSRK3     "3"
809: #define TSRK3BS   "3bs"
810: #define TSRK4     "4"
811: #define TSRK5F    "5f"
812: #define TSRK5DP   "5dp"
813: #define TSRK5BS   "5bs"
814: #define TSRK6VR   "6vr"
815: #define TSRK7VR   "7vr"
816: #define TSRK8VR   "8vr"

818: PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
819: PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
820: PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS,PetscBool);
821: PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS,PetscBool*);
822: PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
823: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
824: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
825: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
826: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);

828: /*J
829:    TSMPRKType - String with the name of a Partitioned Runge-Kutta method

831:    Level: beginner

833: .seealso: TSMPRKSetType(), TS, TSMPRK, TSMPRKRegister()
834: J*/
835: typedef const char* TSMPRKType;
836: #define TSMPRK2A22   "2a22"
837: #define TSMPRK2A23   "2a23"
838: #define TSMPRK2A32   "2a32"
839: #define TSMPRK2A33   "2a33"
840: #define TSMPRKP2    "p2"
841: #define TSMPRKP3    "p3"

843: PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType*);
844: PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType);
845: PETSC_EXTERN PetscErrorCode TSMPRKRegister(TSMPRKType,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscInt[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscInt [],const PetscReal[],const PetscReal[],const PetscReal[]);
846: PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
847: PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
848: PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);

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

853:    Level: beginner

855: .seealso: TSGLEESetType(), TS, TSGLEE, TSGLEERegister()
856: J*/
857: typedef const char* TSGLEEType;
858: #define TSGLEEi1      "BE1"
859: #define TSGLEE23      "23"
860: #define TSGLEE24      "24"
861: #define TSGLEE25I     "25i"
862: #define TSGLEE35      "35"
863: #define TSGLEEEXRK2A  "exrk2a"
864: #define TSGLEERK32G1  "rk32g1"
865: #define TSGLEERK285EX "rk285ex"
866: /*J
867:     TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method.

869:    Level: beginner

871: .seealso: TSGLEESetMode(), TS, TSGLEE, TSGLEERegister()
872: J*/
873: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType*);
874: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts,TSGLEEType);
875: 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[]);
876: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
877: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
878: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);

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

883:    Level: beginner

885: .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
886: J*/
887: typedef const char* TSARKIMEXType;
888: #define TSARKIMEX1BEE   "1bee"
889: #define TSARKIMEXA2     "a2"
890: #define TSARKIMEXL2     "l2"
891: #define TSARKIMEXARS122 "ars122"
892: #define TSARKIMEX2C     "2c"
893: #define TSARKIMEX2D     "2d"
894: #define TSARKIMEX2E     "2e"
895: #define TSARKIMEXPRSSP2 "prssp2"
896: #define TSARKIMEX3      "3"
897: #define TSARKIMEXBPR3   "bpr3"
898: #define TSARKIMEXARS443 "ars443"
899: #define TSARKIMEX4      "4"
900: #define TSARKIMEX5      "5"
901: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
902: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
903: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
904: 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[]);
905: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
906: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
907: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);

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

912:    Level: beginner

914: .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
915: J*/
916: typedef const char* TSRosWType;
917: #define TSROSW2M          "2m"
918: #define TSROSW2P          "2p"
919: #define TSROSWRA3PW       "ra3pw"
920: #define TSROSWRA34PW2     "ra34pw2"
921: #define TSROSWRODAS3      "rodas3"
922: #define TSROSWSANDU3      "sandu3"
923: #define TSROSWASSP3P3S1C  "assp3p3s1c"
924: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
925: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
926: #define TSROSWARK3        "ark3"
927: #define TSROSWTHETA1      "theta1"
928: #define TSROSWTHETA2      "theta2"
929: #define TSROSWGRK4T       "grk4t"
930: #define TSROSWSHAMP4      "shamp4"
931: #define TSROSWVELDD4      "veldd4"
932: #define TSROSW4L          "4l"

934: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS,TSRosWType*);
935: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS,TSRosWType);
936: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
937: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
938: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
939: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
940: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
941: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);

943: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS,PetscInt);
944: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS,PetscInt*);

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

949:   Level: beginner

951:   .seealso: TSBasicSymplecticSetType(), TS, TSBASICSYMPLECTIC, TSBasicSymplecticRegister()
952: J*/
953: typedef const char* TSBasicSymplecticType;
954: #define TSBASICSYMPLECTICSIEULER   "1"
955: #define TSBASICSYMPLECTICVELVERLET "2"
956: #define TSBASICSYMPLECTIC3         "3"
957: #define TSBASICSYMPLECTIC4         "4"
958: PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS,TSBasicSymplecticType);
959: PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS,TSBasicSymplecticType*);
960: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType,PetscInt,PetscInt,PetscReal[],PetscReal[]);
961: PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
962: PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
963: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);

965: /*
966:        PETSc interface to Sundials
967: */
968: #ifdef PETSC_HAVE_SUNDIALS
969: typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
970: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
971: typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
972: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
973: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
974: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
975: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
976: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
977: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
978: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
979: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
980: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
981: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
982: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
983: PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
984: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
985: #endif

987: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
988: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
989: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
990: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);

992: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
993: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
994: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);

996: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS,PetscReal);
997: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS,PetscReal,PetscReal,PetscReal,PetscReal);
998: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*);

1000: PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
1001: PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);

1003: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
1004: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);

1006: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS,PetscBool*);
1007: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS,PetscBool*);

1009: PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
1010: PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
1011: PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
1012: PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1013: PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1014: PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1015: PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst);
1016: #endif