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