Actual source code: petscts.h
petsc-3.12.5 2020-03-29
1: /*
2: User interface for the timestepping package. This package
3: is for use in solving time-dependent PDEs.
4: */
5: #if !defined(PETSCTS_H)
6: #define PETSCTS_H
7: #include <petscsnes.h>
9: /*S
10: TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
12: Level: beginner
14: .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
15: S*/
16: typedef struct _p_TS* TS;
18: /*J
19: TSType - String with the name of a PETSc TS method.
21: Level: beginner
23: .seealso: TSSetType(), TS, TSRegister()
24: J*/
25: typedef const char* TSType;
26: #define TSEULER "euler"
27: #define TSBEULER "beuler"
28: #define TSBASICSYMPLECTIC "basicsymplectic"
29: #define TSPSEUDO "pseudo"
30: #define TSCN "cn"
31: #define TSSUNDIALS "sundials"
32: #define TSRK "rk"
33: #define TSPYTHON "python"
34: #define TSTHETA "theta"
35: #define TSALPHA "alpha"
36: #define TSALPHA2 "alpha2"
37: #define TSGLLE "glle"
38: #define TSGLEE "glee"
39: #define TSSSP "ssp"
40: #define TSARKIMEX "arkimex"
41: #define TSROSW "rosw"
42: #define TSEIMEX "eimex"
43: #define TSMIMEX "mimex"
44: #define TSBDF "bdf"
45: #define TSRADAU5 "radau5"
46: #define TSMPRK "mprk"
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:
63: this must match petsc/finclude/petscts.h
65: Supported types are:
66: TS_EQ_UNSPECIFIED (default)
67: TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
68: TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
70: .seealso: TSGetEquationType(), TSSetEquationType()
71: E*/
72: typedef enum {
73: TS_EQ_UNSPECIFIED = -1,
74: TS_EQ_EXPLICIT = 0,
75: TS_EQ_ODE_EXPLICIT = 1,
76: TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100,
77: TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200,
78: TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300,
79: TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
80: TS_EQ_IMPLICIT = 1000,
81: TS_EQ_ODE_IMPLICIT = 1001,
82: TS_EQ_DAE_IMPLICIT_INDEX1 = 1100,
83: TS_EQ_DAE_IMPLICIT_INDEX2 = 1200,
84: TS_EQ_DAE_IMPLICIT_INDEX3 = 1300,
85: TS_EQ_DAE_IMPLICIT_INDEXHI = 1500
86: } TSEquationType;
87: PETSC_EXTERN const char *const*TSEquationTypes;
89: /*E
90: TSConvergedReason - reason a TS method has converged or not
92: Level: beginner
94: Developer Notes:
95: this must match petsc/finclude/petscts.h
97: Each reason has its own manual page.
99: .seealso: TSGetConvergedReason()
100: E*/
101: typedef enum {
102: TS_CONVERGED_ITERATING = 0,
103: TS_CONVERGED_TIME = 1,
104: TS_CONVERGED_ITS = 2,
105: TS_CONVERGED_USER = 3,
106: TS_CONVERGED_EVENT = 4,
107: TS_CONVERGED_PSEUDO_FATOL = 5,
108: TS_CONVERGED_PSEUDO_FRTOL = 6,
109: TS_DIVERGED_NONLINEAR_SOLVE = -1,
110: TS_DIVERGED_STEP_REJECTED = -2,
111: TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
112: TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
113: } TSConvergedReason;
114: PETSC_EXTERN const char *const*TSConvergedReasons;
116: /*MC
117: TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
119: Level: beginner
121: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
122: M*/
124: /*MC
125: TS_CONVERGED_TIME - the final time was reached
127: Level: beginner
129: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxTime(), TSGetMaxTime(), TSGetSolveTime()
130: M*/
132: /*MC
133: TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
135: Level: beginner
137: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxSteps(), TSGetMaxSteps()
138: M*/
140: /*MC
141: TS_CONVERGED_USER - user requested termination
143: Level: beginner
145: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
146: M*/
148: /*MC
149: TS_CONVERGED_EVENT - user requested termination on event detection
151: Level: beginner
153: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
154: M*/
156: /*MC
157: TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for TSPSEUDO
159: Level: beginner
161: Options Database:
162: . -ts_pseudo_frtol <rtol>
164: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FATOL
165: M*/
167: /*MC
168: TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for TSPSEUDO
170: Level: beginner
172: Options Database:
173: . -ts_pseudo_fatol <atol>
175: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FRTOL
176: M*/
178: /*MC
179: TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
181: Level: beginner
183: Notes:
184: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
186: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
187: M*/
189: /*MC
190: TS_DIVERGED_STEP_REJECTED - too many steps were rejected
192: Level: beginner
194: Notes:
195: See TSSetMaxStepRejections() for how to allow more step rejections.
197: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
198: M*/
200: /*E
201: TSExactFinalTimeOption - option for handling of final time step
203: Level: beginner
205: Developer Notes:
206: this must match petsc/finclude/petscts.h
208: $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded
209: $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
210: $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
212: .seealso: TSGetConvergedReason(), TSSetExactFinalTime(), TSGetExactFinalTime()
214: E*/
215: typedef enum {TS_EXACTFINALTIME_UNSPECIFIED=0,TS_EXACTFINALTIME_STEPOVER=1,TS_EXACTFINALTIME_INTERPOLATE=2,TS_EXACTFINALTIME_MATCHSTEP=3} TSExactFinalTimeOption;
216: PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
218: /* Logging support */
219: PETSC_EXTERN PetscClassId TS_CLASSID;
220: PETSC_EXTERN PetscClassId DMTS_CLASSID;
221: PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
223: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
224: PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
226: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
227: PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
228: PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
230: PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
231: PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
232: PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
233: PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
234: PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
235: PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
237: PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
238: PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
239: PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
240: PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
241: PETSC_EXTERN PetscErrorCode TSSetUp(TS);
242: PETSC_EXTERN PetscErrorCode TSReset(TS);
244: PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
245: PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
247: PETSC_EXTERN PetscErrorCode TS2SetSolution(TS,Vec,Vec);
248: PETSC_EXTERN PetscErrorCode TS2GetSolution(TS,Vec*,Vec*);
250: PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS,PetscInt*,Vec*);
251: PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS,Vec*);
252: PETSC_EXTERN PetscErrorCode TSGetTimeError(TS,PetscInt,Vec*);
253: PETSC_EXTERN PetscErrorCode TSSetTimeError(TS,Vec);
255: PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
256: PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS,Mat*,PetscErrorCode(**)(TS,PetscReal,Vec,Mat,void*),void**);
257: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS,PetscReal,Vec,Mat);
258: PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void*);
259: PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS,PetscReal,Vec,Vec,PetscReal,Mat,PetscBool);
260: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
261: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSComputeDRDUFunction(TS,PetscReal,Vec,Vec*);
262: PETSC_EXTERN PetscErrorCode TSSetIHessianProduct(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
263: 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: void*);
267: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
268: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
269: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
270: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
271: PETSC_EXTERN PetscErrorCode TSSetRHSHessianProduct(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
272: 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: void*);
276: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
277: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
278: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
279: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
280: PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS,PetscInt,Vec*,Vec*,Vec);
281: PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS,PetscInt*,Vec**,Vec**,Vec*);
282: PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS,Vec,Mat,Mat);
284: /*S
285: TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
287: Level: advanced
289: .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectorySetType(), TSTrajectoryDestroy(), TSTrajectoryReset()
290: S*/
291: typedef struct _p_TSTrajectory* TSTrajectory;
293: /*J
294: TSTrajectorySetType - String with the name of a PETSc TS trajectory storage method
296: Level: intermediate
298: .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy()
299: J*/
300: typedef const char* TSTrajectoryType;
301: #define TSTRAJECTORYBASIC "basic"
302: #define TSTRAJECTORYSINGLEFILE "singlefile"
303: #define TSTRAJECTORYMEMORY "memory"
304: #define TSTRAJECTORYVISUALIZATION "visualization"
306: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
307: PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID;
308: PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled;
310: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
311: PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
313: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
314: PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
315: PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
316: PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory,PetscViewer);
317: PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,TSTrajectoryType);
318: PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory,TS,TSTrajectoryType*);
319: PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
320: PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
321: PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory,TS,PetscInt,PetscReal*,Vec,Vec);
322: PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory,TS,PetscReal,Vec*,Vec*);
323: PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory,PetscInt*);
324: PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory,Vec*,Vec*);
325: PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
326: PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[],PetscErrorCode (*)(TSTrajectory,TS));
327: PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
328: PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
329: PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory,PetscBool);
330: PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory,PetscBool);
331: PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory,const char * const*);
332: PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
333: PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory,PetscBool);
334: PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory,PetscBool*);
335: PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory,PetscBool);
336: PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory,const char[]);
337: PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory,const char[]);
338: PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS,TSTrajectory*);
340: PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
341: PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
342: 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*);
343: PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
344: PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec);
345: PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS,PetscBool,TS*);
346: PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS,PetscBool*,TS*);
348: PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems*,TS);
349: PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
350: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
351: PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
352: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
354: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetRHSJacobianP()") PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
355: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSComputeRHSJacobianP()") PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
356: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
357: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
358: PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
359: PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);
361: PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
362: PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
363: PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
364: PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
365: PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS,Mat);
366: PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
368: PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS,PetscInt,Mat);
369: PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS,PetscInt*,Mat*);
370: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSCreateQuadratureTS() and TSForwardSetSensitivities() (since version 3.12)") PetscErrorCode TSForwardSetIntegralGradients(TS,PetscInt,Vec *);
371: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSForwardGetSensitivities()") PetscErrorCode TSForwardGetIntegralGradients(TS,PetscInt*,Vec **);
372: PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
373: PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
374: PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
375: PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
376: PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS,Mat);
377: PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS,PetscInt*,Mat**);
379: PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS,PetscInt);
380: PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS,PetscInt*);
381: PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS,PetscReal);
382: PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS,PetscReal*);
383: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
384: PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS,TSExactFinalTimeOption*);
386: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetTime[Step]() (since version 3.8)") PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
387: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
388: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
389: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)") PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
390: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)") PetscErrorCode TSGetTotalSteps(TS,PetscInt*);
392: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
393: PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
395: typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx;
396: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
397: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
398: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
399: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
400: PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
401: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS,PetscInt,PetscReal,Vec,void*);
403: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat *);
404: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
406: PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
407: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
408: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
410: PETSC_EXTERN PetscErrorCode TSStep(TS);
411: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*);
412: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
413: PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
414: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
415: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
416: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
417: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
418: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
419: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
420: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
421: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
422: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
423: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
424: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
425: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
426: PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
427: PETSC_EXTERN PetscErrorCode TSRollBack(TS);
429: PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
431: PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
432: PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
433: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
434: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
435: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
436: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS,PetscInt*);
437: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS,PetscInt);
439: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
440: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
441: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS,PetscReal,Vec,Mat,void*);
442: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
443: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
444: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
445: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
446: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
448: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
449: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
450: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*);
451: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,TSForcingFunction,void*);
453: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
454: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
455: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
456: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
457: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
458: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
460: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS,PetscReal,Vec,Vec,Vec,Vec,void*);
461: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat,void*);
462: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS,Vec,TSI2Function,void*);
463: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS,Vec*,TSI2Function*,void**);
464: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS,Mat,Mat,TSI2Jacobian,void*);
465: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS,Mat*,Mat*,TSI2Jacobian*,void**);
467: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS,const char[],IS);
468: PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS,const char[],IS*);
469: PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS,const char[],Vec,TSRHSFunction,void*);
470: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS,const char[],TS*);
471: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS,PetscInt*,TS*[]);
472: PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
473: PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool*);
475: PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
476: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
477: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
478: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
479: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
480: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
481: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
483: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
484: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
485: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
486: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode(*)(TS));
487: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
488: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
489: PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
490: PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
491: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
492: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
493: PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
494: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
495: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
496: PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
497: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
498: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
499: PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
500: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
501: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS,Vec,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
502: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
503: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
504: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
505: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);
507: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
508: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
509: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
510: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
511: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
512: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
513: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
514: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
515: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
517: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
519: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
520: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
521: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
522: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
523: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS,PetscReal,Vec,Vec,Vec,Vec);
524: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat);
525: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
527: PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
529: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
530: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
531: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
532: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
533: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
534: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
535: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
536: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
537: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
538: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM,TSI2Function,void*);
539: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM,TSI2Function*,void**);
540: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM,TSI2Jacobian,void*);
541: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM,TSI2Jacobian*,void**);
543: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
544: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
545: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,TSForcingFunction,void*);
546: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,TSForcingFunction*,void**);
547: PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
548: PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
549: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **);
551: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
552: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
553: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
555: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
556: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
558: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
559: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
560: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
561: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
563: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
564: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
565: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
566: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
568: PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
569: PETSC_EXTERN PetscErrorCode DMPlexTSGetGradientDM(DM,PetscFV,DM*);
571: typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx;
572: typedef struct {
573: Vec ray;
574: VecScatter scatter;
575: PetscViewer viewer;
576: TSMonitorLGCtx lgctx;
577: } TSMonitorDMDARayCtx;
578: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
579: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
580: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
582: /* Dynamic creation and loading functions */
583: PETSC_EXTERN PetscFunctionList TSList;
584: PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
585: PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
586: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
588: PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
589: PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
590: PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
592: PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
593: PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
594: PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
595: PETSC_STATIC_INLINE PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
597: #define TS_FILE_CLASSID 1211225
599: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
600: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
602: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
603: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
604: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
605: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
606: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
607: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
608: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
609: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
610: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
611: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
612: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
613: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
614: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
615: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
616: PETSC_EXTERN PetscErrorCode TSMonitorError(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat *);
618: typedef struct _n_TSMonitorEnvelopeCtx* TSMonitorEnvelopeCtx;
619: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
620: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
621: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
622: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);
624: typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx;
625: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
626: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
627: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
629: typedef struct _n_TSMonitorSPCtx* TSMonitorSPCtx;
630: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPCtx*);
631: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx*);
632: PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS,PetscInt,PetscReal,Vec,void*);
634: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
635: PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS,PetscReal);
636: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]);
637: /*J
638: TSSSPType - string with the name of TSSSP scheme.
640: Level: beginner
642: .seealso: TSSSPSetType(), TS
643: J*/
644: typedef const char* TSSSPType;
645: #define TSSSPRKS2 "rks2"
646: #define TSSSPRKS3 "rks3"
647: #define TSSSPRK104 "rk104"
649: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
650: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
651: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
652: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
653: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
654: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
655: PETSC_EXTERN PetscFunctionList TSSSPList;
657: /*S
658: TSAdapt - Abstract object that manages time-step adaptivity
660: Level: beginner
662: .seealso: TS, TSAdaptCreate(), TSAdaptType
663: S*/
664: typedef struct _p_TSAdapt *TSAdapt;
666: /*E
667: TSAdaptType - String with the name of TSAdapt scheme.
669: Level: beginner
671: .seealso: TSAdaptSetType(), TS
672: E*/
673: typedef const char *TSAdaptType;
674: #define TSADAPTNONE "none"
675: #define TSADAPTBASIC "basic"
676: #define TSADAPTDSP "dsp"
677: #define TSADAPTCFL "cfl"
678: #define TSADAPTGLEE "glee"
679: #define TSADAPTHISTORY "history"
681: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
682: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
683: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
684: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
685: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
686: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
687: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt,TSAdaptType*);
688: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
689: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
690: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
691: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
692: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
693: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
694: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
695: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
696: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
697: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
698: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
699: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
700: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt,PetscBool);
701: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt,PetscReal,PetscReal);
702: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt,PetscReal*,PetscReal*);
703: PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt,PetscReal);
704: PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt,PetscReal*);
705: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt,PetscReal,PetscReal);
706: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt,PetscReal*,PetscReal*);
707: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
708: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt,PetscReal*,PetscReal*);
709: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
710: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt,PetscInt n,PetscReal hist[],PetscBool);
711: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt,TSTrajectory,PetscBool);
712: PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt,PetscInt,PetscReal*,PetscReal*);
713: PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt,PetscInt);
714: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt,const char *);
715: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt,PetscReal,PetscReal,PetscReal);
717: /*S
718: TSGLLEAdapt - Abstract object that manages time-step adaptivity
720: Level: beginner
722: Developer Notes:
723: This functionality should be replaced by the TSAdapt.
725: .seealso: TSGLLE, TSGLLEAdaptCreate(), TSGLLEAdaptType
726: S*/
727: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
729: /*J
730: TSGLLEAdaptType - String with the name of TSGLLEAdapt scheme
732: Level: beginner
734: .seealso: TSGLLEAdaptSetType(), TS
735: J*/
736: typedef const char *TSGLLEAdaptType;
737: #define TSGLLEADAPT_NONE "none"
738: #define TSGLLEADAPT_SIZE "size"
739: #define TSGLLEADAPT_BOTH "both"
741: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[],PetscErrorCode (*)(TSGLLEAdapt));
742: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
743: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
744: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm,TSGLLEAdapt*);
745: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt,TSGLLEAdaptType);
746: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt,const char[]);
747: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
748: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt,PetscViewer);
749: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems*,TSGLLEAdapt);
750: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt*);
752: /*J
753: TSGLLEAcceptType - String with the name of TSGLLEAccept scheme
755: Level: beginner
757: .seealso: TSGLLESetAcceptType(), TS
758: J*/
759: typedef const char *TSGLLEAcceptType;
760: #define TSGLLEACCEPT_ALWAYS "always"
762: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
763: PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[],TSGLLEAcceptFunction);
765: /*J
766: TSGLLEType - family of time integration method within the General Linear class
768: Level: beginner
770: .seealso: TSGLLESetType(), TSGLLERegister()
771: J*/
772: typedef const char* TSGLLEType;
773: #define TSGLLE_IRKS "irks"
775: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[],PetscErrorCode(*)(TS));
776: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
777: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
778: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS,TSGLLEType);
779: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS,TSGLLEAdapt*);
780: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS,TSGLLEAcceptType);
782: /*J
783: TSEIMEXType - String with the name of an Extrapolated IMEX method.
785: Level: beginner
787: .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
788: J*/
789: #define TSEIMEXType char*
791: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
792: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
793: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
795: /*J
796: TSRKType - String with the name of a Runge-Kutta method.
798: Level: beginner
800: .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
801: J*/
802: typedef const char* TSRKType;
803: #define TSRK1FE "1fe"
804: #define TSRK2A "2a"
805: #define TSRK3 "3"
806: #define TSRK3BS "3bs"
807: #define TSRK4 "4"
808: #define TSRK5F "5f"
809: #define TSRK5DP "5dp"
810: #define TSRK5BS "5bs"
811: #define TSRK6VR "6vr"
812: #define TSRK7VR "7vr"
813: #define TSRK8VR "8vr"
815: PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
816: PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
817: PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS,PetscBool);
818: PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS,PetscBool*);
819: PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
820: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
821: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
822: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
823: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
825: /*J
826: TSMPRKType - String with the name of a Partitioned Runge-Kutta method
828: Level: beginner
830: .seealso: TSMPRKSetType(), TS, TSMPRK, TSMPRKRegister()
831: J*/
832: typedef const char* TSMPRKType;
833: #define TSMPRK2A22 "2a22"
834: #define TSMPRK2A23 "2a23"
835: #define TSMPRK2A32 "2a32"
836: #define TSMPRK2A33 "2a33"
837: #define TSMPRKP2 "p2"
838: #define TSMPRKP3 "p3"
840: PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType*);
841: PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType);
842: 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[]);
843: PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
844: PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
845: PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
847: /*J
848: TSGLEEType - String with the name of a General Linear with Error Estimation method.
850: Level: beginner
852: .seealso: TSGLEESetType(), TS, TSGLEE, TSGLEERegister()
853: J*/
854: typedef const char* TSGLEEType;
855: #define TSGLEEi1 "BE1"
856: #define TSGLEE23 "23"
857: #define TSGLEE24 "24"
858: #define TSGLEE25I "25i"
859: #define TSGLEE35 "35"
860: #define TSGLEEEXRK2A "exrk2a"
861: #define TSGLEERK32G1 "rk32g1"
862: #define TSGLEERK285EX "rk285ex"
863: /*J
864: TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method.
866: Level: beginner
868: .seealso: TSGLEESetMode(), TS, TSGLEE, TSGLEERegister()
869: J*/
870: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType*);
871: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts,TSGLEEType);
872: 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[]);
873: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
874: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
875: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
877: /*J
878: TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
880: Level: beginner
882: .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
883: J*/
884: typedef const char* TSARKIMEXType;
885: #define TSARKIMEX1BEE "1bee"
886: #define TSARKIMEXA2 "a2"
887: #define TSARKIMEXL2 "l2"
888: #define TSARKIMEXARS122 "ars122"
889: #define TSARKIMEX2C "2c"
890: #define TSARKIMEX2D "2d"
891: #define TSARKIMEX2E "2e"
892: #define TSARKIMEXPRSSP2 "prssp2"
893: #define TSARKIMEX3 "3"
894: #define TSARKIMEXBPR3 "bpr3"
895: #define TSARKIMEXARS443 "ars443"
896: #define TSARKIMEX4 "4"
897: #define TSARKIMEX5 "5"
898: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
899: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
900: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
901: 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[]);
902: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
903: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
904: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
906: /*J
907: TSRosWType - String with the name of a Rosenbrock-W method.
909: Level: beginner
911: .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
912: J*/
913: typedef const char* TSRosWType;
914: #define TSROSW2M "2m"
915: #define TSROSW2P "2p"
916: #define TSROSWRA3PW "ra3pw"
917: #define TSROSWRA34PW2 "ra34pw2"
918: #define TSROSWRODAS3 "rodas3"
919: #define TSROSWSANDU3 "sandu3"
920: #define TSROSWASSP3P3S1C "assp3p3s1c"
921: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
922: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
923: #define TSROSWARK3 "ark3"
924: #define TSROSWTHETA1 "theta1"
925: #define TSROSWTHETA2 "theta2"
926: #define TSROSWGRK4T "grk4t"
927: #define TSROSWSHAMP4 "shamp4"
928: #define TSROSWVELDD4 "veldd4"
929: #define TSROSW4L "4l"
931: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS,TSRosWType*);
932: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS,TSRosWType);
933: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
934: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
935: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
936: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
937: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
938: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
940: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS,PetscInt);
941: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS,PetscInt*);
943: /*J
944: TSBasicSymplecticType - String with the name of a basic symplectic integration method.
946: Level: beginner
948: .seealso: TSBasicSymplecticSetType(), TS, TSBASICSYMPLECTIC, TSBasicSymplecticRegister()
949: J*/
950: typedef const char* TSBasicSymplecticType;
951: #define TSBASICSYMPLECTICSIEULER "1"
952: #define TSBASICSYMPLECTICVELVERLET "2"
953: #define TSBASICSYMPLECTIC3 "3"
954: #define TSBASICSYMPLECTIC4 "4"
955: PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS,TSBasicSymplecticType);
956: PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS,TSBasicSymplecticType*);
957: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType,PetscInt,PetscInt,PetscReal[],PetscReal[]);
958: PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
959: PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
960: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
962: /*
963: PETSc interface to Sundials
964: */
965: #ifdef PETSC_HAVE_SUNDIALS
966: typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
967: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
968: typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
969: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
970: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
971: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
972: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
973: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
974: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
975: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
976: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
977: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
978: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
979: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
980: PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
981: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
982: #endif
984: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
985: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
986: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
987: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
989: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
990: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
991: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
993: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS,PetscReal);
994: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS,PetscReal,PetscReal,PetscReal,PetscReal);
995: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
997: PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
998: PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
1000: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
1001: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
1003: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS,PetscBool*);
1004: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS,PetscBool*);
1005: #endif