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