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