Actual source code: petscts.h
petsc-3.8.0 2017-09-26
1: /*
2: User interface for the timestepping package. This package
3: is for use in solving time-dependent PDEs.
4: */
7: #include <petscsnes.h>
9: /*S
10: TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
12: Level: beginner
14: Concepts: ODE solvers
16: .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
17: S*/
18: typedef struct _p_TS* TS;
20: /*J
21: TSType - String with the name of a PETSc TS method.
23: Level: beginner
25: .seealso: TSSetType(), TS, TSRegister()
26: J*/
27: typedef const char* TSType;
28: #define TSEULER "euler"
29: #define TSBEULER "beuler"
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"
47: /*E
48: TSProblemType - Determines the type of problem this TS object is to be used to solve
50: Level: beginner
52: .seealso: TSCreate()
53: E*/
54: typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
56: /*E
57: TSEquationType - type of TS problem that is solved
59: Level: beginner
61: Developer Notes: this must match petsc/finclude/petscts.h
63: Supported types are:
64: TS_EQ_UNSPECIFIED (default)
65: TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
66: TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
68: .seealso: TSGetEquationType(), TSSetEquationType()
69: E*/
70: typedef enum {
71: TS_EQ_UNSPECIFIED = -1,
72: TS_EQ_EXPLICIT = 0,
73: TS_EQ_ODE_EXPLICIT = 1,
74: TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100,
75: TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200,
76: TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300,
77: TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
78: TS_EQ_IMPLICIT = 1000,
79: TS_EQ_ODE_IMPLICIT = 1001,
80: TS_EQ_DAE_IMPLICIT_INDEX1 = 1100,
81: TS_EQ_DAE_IMPLICIT_INDEX2 = 1200,
82: TS_EQ_DAE_IMPLICIT_INDEX3 = 1300,
83: TS_EQ_DAE_IMPLICIT_INDEXHI = 1500
84: } TSEquationType;
85: PETSC_EXTERN const char *const*TSEquationTypes;
87: /*E
88: TSConvergedReason - reason a TS method has converged or not
90: Level: beginner
92: Developer Notes: this must match petsc/finclude/petscts.h
94: Each reason has its own manual page.
96: .seealso: TSGetConvergedReason()
97: E*/
98: typedef enum {
99: TS_CONVERGED_ITERATING = 0,
100: TS_CONVERGED_TIME = 1,
101: TS_CONVERGED_ITS = 2,
102: TS_CONVERGED_USER = 3,
103: TS_CONVERGED_EVENT = 4,
104: TS_CONVERGED_PSEUDO_FATOL = 5,
105: TS_CONVERGED_PSEUDO_FRTOL = 6,
106: TS_DIVERGED_NONLINEAR_SOLVE = -1,
107: TS_DIVERGED_STEP_REJECTED = -2
108: } TSConvergedReason;
109: PETSC_EXTERN const char *const*TSConvergedReasons;
111: /*MC
112: TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
114: Level: beginner
116: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
117: M*/
119: /*MC
120: TS_CONVERGED_TIME - the final time was reached
122: Level: beginner
124: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxTime(), TSGetMaxTime(), TSGetSolveTime()
125: M*/
127: /*MC
128: TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
130: Level: beginner
132: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxSteps(), TSGetMaxSteps()
133: M*/
135: /*MC
136: TS_CONVERGED_USER - user requested termination
138: Level: beginner
140: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
141: M*/
143: /*MC
144: TS_CONVERGED_EVENT - user requested termination on event detection
146: Level: beginner
148: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
149: M*/
151: /*MC
152: TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for TSPSEUDO
154: Level: beginner
156: Options Database:
157: . -ts_pseudo_frtol <rtol>
159: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FATOL
160: M*/
162: /*MC
163: TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for TSPSEUDO
165: Level: beginner
167: Options Database:
168: . -ts_pseudo_fatol <atol>
170: .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FRTOL
171: M*/
173: /*MC
174: TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
176: Level: beginner
178: Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
180: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
181: M*/
183: /*MC
184: TS_DIVERGED_STEP_REJECTED - too many steps were rejected
186: Level: beginner
188: Notes: See TSSetMaxStepRejections() for how to allow more step rejections.
190: .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
191: M*/
193: /*E
194: TSExactFinalTimeOption - option for handling of final time step
196: Level: beginner
198: Developer Notes: this must match petsc/finclude/petscts.h
200: $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded
201: $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
202: $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
204: .seealso: TSGetConvergedReason(), TSSetExactFinalTime(), TSGetExactFinalTime()
206: E*/
207: typedef enum {TS_EXACTFINALTIME_UNSPECIFIED=0,TS_EXACTFINALTIME_STEPOVER=1,TS_EXACTFINALTIME_INTERPOLATE=2,TS_EXACTFINALTIME_MATCHSTEP=3} TSExactFinalTimeOption;
208: PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
210: /* Logging support */
211: PETSC_EXTERN PetscClassId TS_CLASSID;
212: PETSC_EXTERN PetscClassId DMTS_CLASSID;
213: PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
215: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
216: PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
218: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
219: PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
220: PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
222: PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
223: PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
224: PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
225: PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
226: PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
227: PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
229: PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
230: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
231: PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
232: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
234: PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
235: PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
236: PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
237: PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
238: PETSC_EXTERN PetscErrorCode TSSetUp(TS);
239: PETSC_EXTERN PetscErrorCode TSReset(TS);
241: PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
242: PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
244: PETSC_EXTERN PetscErrorCode TS2SetSolution(TS,Vec,Vec);
245: PETSC_EXTERN PetscErrorCode TS2GetSolution(TS,Vec*,Vec*);
247: PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS,PetscInt*,Vec*);
248: PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS,Vec*);
249: PETSC_EXTERN PetscErrorCode TSGetTimeError(TS,PetscInt,Vec*);
250: PETSC_EXTERN PetscErrorCode TSSetTimeError(TS,Vec);
252: /*S
253: TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/ADE at each time step)
255: Level: advanced
257: Concepts: ODE solvers, trajectory
259: .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectorySetType(), TSTrajectoryDestroy()
260: S*/
261: typedef struct _p_TSTrajectory* TSTrajectory;
263: /*J
264: TSTrajectorySetType - String with the name of a PETSc TS trajectory storage method
266: Level: intermediate
268: .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy()
269: J*/
270: typedef const char* TSTrajectoryType;
271: #define TSTRAJECTORYBASIC "basic"
272: #define TSTRAJECTORYSINGLEFILE "singlefile"
273: #define TSTRAJECTORYMEMORY "memory"
274: #define TSTRAJECTORYVISUALIZATION "visualization"
276: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
277: PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID;
278: PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled;
280: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
282: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
283: PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
284: PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory,PetscViewer);
285: PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,const TSTrajectoryType);
286: PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
287: PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
288: PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
289: PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[],PetscErrorCode (*)(TSTrajectory,TS));
290: PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
291: PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
292: PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory,PetscBool);
293: PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory,const char * const*);
294: PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
295: PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS,TSTrajectory*);
297: PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
298: PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
299: PETSC_EXTERN 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*);
300: PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
301: PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec);
303: PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
304: PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
305: PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);
307: PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
308: PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
309: PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
310: PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
311: PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
312: PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
314: PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS,PetscInt,Vec*,PetscInt,Vec*);
315: PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS,PetscInt*,Vec**,PetscInt*,Vec**);
316: PETSC_EXTERN PetscErrorCode TSForwardSetIntegralGradients(TS,PetscInt,Vec *,Vec *);
317: PETSC_EXTERN PetscErrorCode TSForwardGetIntegralGradients(TS,PetscInt*,Vec **,Vec **);
318: PETSC_EXTERN PetscErrorCode TSForwardSetRHSJacobianP(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,void*),void*);
319: PETSC_EXTERN PetscErrorCode TSForwardComputeRHSJacobianP(TS,PetscReal,Vec,Vec*);
320: PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
321: PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
322: PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
324: PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS,PetscInt);
325: PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS,PetscInt*);
326: PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS,PetscReal);
327: PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS,PetscReal*);
328: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
329: PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS,TSExactFinalTimeOption*);
331: PETSC_EXTERN PETSC_DEPRECATED("Use TSSetTime[Step]") PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
332: PETSC_EXTERN PETSC_DEPRECATED("Use TSSetMax{Steps|Time}") PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
333: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetMax{Steps|Time}") PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
334: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber") PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
335: PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber") PetscErrorCode TSGetTotalSteps(TS,PetscInt*);
337: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
339: typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx;
340: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
341: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
342: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
343: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
344: PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
346: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat *);
347: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
349: PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
350: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
351: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
353: PETSC_EXTERN PetscErrorCode TSStep(TS);
354: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*);
355: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
356: PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
357: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
358: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
359: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
360: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
361: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
362: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
363: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
364: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
365: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
366: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
367: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
368: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
369: PETSC_EXTERN PetscErrorCode TSRollBack(TS);
371: PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
373: PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
374: PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
375: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
376: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
377: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
378: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS,PetscInt*);
379: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS,PetscInt);
381: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
382: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
383: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
384: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
385: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
386: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
387: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
389: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
390: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
391: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*);
392: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,TSForcingFunction,void*);
394: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
395: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
396: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
397: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
398: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
399: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
401: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS,PetscReal,Vec,Vec,Vec,Vec,void*);
402: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat,void*);
403: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS,Vec,TSI2Function,void*);
404: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS,Vec*,TSI2Function*,void**);
405: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS,Mat,Mat,TSI2Jacobian,void*);
406: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS,Mat*,Mat*,TSI2Jacobian*,void**);
408: PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
409: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
410: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
411: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
412: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
413: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
414: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
416: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
417: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
418: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
419: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode(*)(TS));
420: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
421: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
422: PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
423: PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
424: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
425: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
426: PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
427: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
428: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
429: PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
430: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
431: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
432: PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
433: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
434: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS,Vec,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
435: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
436: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
437: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
438: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);
440: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
441: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
442: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
443: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
444: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
445: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
446: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
447: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
448: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
450: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
452: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
453: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
454: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
455: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
456: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS,PetscReal,Vec,Vec,Vec,Vec);
457: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat);
458: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
460: PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
462: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
463: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
464: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
465: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
466: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
467: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
468: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
469: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
470: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
471: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM,TSI2Function,void*);
472: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM,TSI2Function*,void**);
473: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM,TSI2Jacobian,void*);
474: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM,TSI2Jacobian*,void**);
476: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
477: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
478: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,TSForcingFunction,void*);
479: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,TSForcingFunction*,void**);
480: PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
481: PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
482: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **);
484: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
485: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
486: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
488: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
489: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
491: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
492: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
493: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
494: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
496: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
497: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
498: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
499: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
501: PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
502: PETSC_EXTERN PetscErrorCode DMPlexTSGetGradientDM(DM,PetscFV,DM*);
504: typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx;
505: typedef struct {
506: Vec ray;
507: VecScatter scatter;
508: PetscViewer viewer;
509: TSMonitorLGCtx lgctx;
510: } TSMonitorDMDARayCtx;
511: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
512: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
513: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
515: /* Dynamic creation and loading functions */
516: PETSC_EXTERN PetscFunctionList TSList;
517: PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
518: PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
519: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
521: PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
522: PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
523: PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
525: PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
526: PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
527: PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
528: PETSC_STATIC_INLINE PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
530: #define TS_FILE_CLASSID 1211225
532: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
533: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
535: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
536: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
537: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
538: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
539: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
540: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
541: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
542: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
543: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
544: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
545: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
546: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
547: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
548: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
550: typedef struct _n_TSMonitorEnvelopeCtx* TSMonitorEnvelopeCtx;
551: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
552: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
553: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
554: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);
556: typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx;
557: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
558: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
559: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
561: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
562: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]);
563: /*J
564: TSSSPType - string with the name of TSSSP scheme.
566: Level: beginner
568: .seealso: TSSSPSetType(), TS
569: J*/
570: typedef const char* TSSSPType;
571: #define TSSSPRKS2 "rks2"
572: #define TSSSPRKS3 "rks3"
573: #define TSSSPRK104 "rk104"
575: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
576: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
577: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
578: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
579: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
580: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
581: PETSC_EXTERN PetscFunctionList TSSSPList;
583: /*S
584: TSAdapt - Abstract object that manages time-step adaptivity
586: Level: beginner
588: .seealso: TS, TSAdaptCreate(), TSAdaptType
589: S*/
590: typedef struct _p_TSAdapt *TSAdapt;
592: /*E
593: TSAdaptType - String with the name of TSAdapt scheme.
595: Level: beginner
597: .seealso: TSAdaptSetType(), TS
598: E*/
599: typedef const char *TSAdaptType;
600: #define TSADAPTNONE "none"
601: #define TSADAPTBASIC "basic"
602: #define TSADAPTCFL "cfl"
603: #define TSADAPTGLEE "glee"
605: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
606: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
607: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
608: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
609: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
610: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
611: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt,TSAdaptType*);
612: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
613: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
614: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
615: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
616: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
617: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
618: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
619: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
620: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
621: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
622: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
623: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
624: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt,PetscBool);
625: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt,PetscReal,PetscReal);
626: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt,PetscReal*,PetscReal*);
627: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt,PetscReal,PetscReal);
628: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt,PetscReal*,PetscReal*);
629: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
630: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt,PetscReal*,PetscReal*);
631: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
633: /*S
634: TSGLLEAdapt - Abstract object that manages time-step adaptivity
636: Level: beginner
638: Developer Notes:
639: This functionality should be replaced by the TSAdapt.
641: .seealso: TSGLLE, TSGLLEAdaptCreate(), TSGLLEAdaptType
642: S*/
643: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
645: /*J
646: TSGLLEAdaptType - String with the name of TSGLLEAdapt scheme
648: Level: beginner
650: .seealso: TSGLLEAdaptSetType(), TS
651: J*/
652: typedef const char *TSGLLEAdaptType;
653: #define TSGLLEADAPT_NONE "none"
654: #define TSGLLEADAPT_SIZE "size"
655: #define TSGLLEADAPT_BOTH "both"
657: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[],PetscErrorCode (*)(TSGLLEAdapt));
658: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
659: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
660: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm,TSGLLEAdapt*);
661: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt,TSGLLEAdaptType);
662: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt,const char[]);
663: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
664: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt,PetscViewer);
665: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems*,TSGLLEAdapt);
666: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt*);
668: /*J
669: TSGLLEAcceptType - String with the name of TSGLLEAccept scheme
671: Level: beginner
673: .seealso: TSGLLESetAcceptType(), TS
674: J*/
675: typedef const char *TSGLLEAcceptType;
676: #define TSGLLEACCEPT_ALWAYS "always"
678: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
679: PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[],TSGLLEAcceptFunction);
681: /*J
682: TSGLLEType - family of time integration method within the General Linear class
684: Level: beginner
686: .seealso: TSGLLESetType(), TSGLLERegister()
687: J*/
688: typedef const char* TSGLLEType;
689: #define TSGLLE_IRKS "irks"
691: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[],PetscErrorCode(*)(TS));
692: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
693: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
694: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS,TSGLLEType);
695: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS,TSGLLEAdapt*);
696: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS,TSGLLEAcceptType);
698: /*J
699: TSEIMEXType - String with the name of an Extrapolated IMEX method.
701: Level: beginner
703: .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
704: J*/
705: #define TSEIMEXType char*
707: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
708: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
709: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
711: /*J
712: TSRKType - String with the name of a Runge-Kutta method.
714: Level: beginner
716: .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
717: J*/
718: typedef const char* TSRKType;
719: #define TSRK1FE "1fe"
720: #define TSRK2A "2a"
721: #define TSRK3 "3"
722: #define TSRK3BS "3bs"
723: #define TSRK4 "4"
724: #define TSRK5F "5f"
725: #define TSRK5DP "5dp"
726: #define TSRK5BS "5bs"
728: PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
729: PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
730: PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
731: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
732: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
733: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
734: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
736: /*J
737: TSGLEEType - String with the name of a General Linear with Error Estimation method.
739: Level: beginner
741: .seealso: TSGLEESetType(), TS, TSGLEE, TSGLEERegister()
742: J*/
743: typedef const char* TSGLEEType;
744: #define TSGLEEi1 "BE1"
745: #define TSGLEE23 "23"
746: #define TSGLEE24 "24"
747: #define TSGLEE25I "25i"
748: #define TSGLEE35 "35"
749: #define TSGLEEEXRK2A "exrk2a"
750: #define TSGLEERK32G1 "rk32g1"
751: #define TSGLEERK285EX "rk285ex"
752: /*J
753: TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method.
755: Level: beginner
757: .seealso: TSGLEESetMode(), TS, TSGLEE, TSGLEERegister()
758: J*/
759: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType*);
760: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts,TSGLEEType);
761: 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[]);
762: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
763: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
764: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
766: /*J
767: TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
769: Level: beginner
771: .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
772: J*/
773: typedef const char* TSARKIMEXType;
774: #define TSARKIMEX1BEE "1bee"
775: #define TSARKIMEXA2 "a2"
776: #define TSARKIMEXL2 "l2"
777: #define TSARKIMEXARS122 "ars122"
778: #define TSARKIMEX2C "2c"
779: #define TSARKIMEX2D "2d"
780: #define TSARKIMEX2E "2e"
781: #define TSARKIMEXPRSSP2 "prssp2"
782: #define TSARKIMEX3 "3"
783: #define TSARKIMEXBPR3 "bpr3"
784: #define TSARKIMEXARS443 "ars443"
785: #define TSARKIMEX4 "4"
786: #define TSARKIMEX5 "5"
787: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
788: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
789: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
790: 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[]);
791: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
792: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
793: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
795: /*J
796: TSRosWType - String with the name of a Rosenbrock-W method.
798: Level: beginner
800: .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
801: J*/
802: typedef const char* TSRosWType;
803: #define TSROSW2M "2m"
804: #define TSROSW2P "2p"
805: #define TSROSWRA3PW "ra3pw"
806: #define TSROSWRA34PW2 "ra34pw2"
807: #define TSROSWRODAS3 "rodas3"
808: #define TSROSWSANDU3 "sandu3"
809: #define TSROSWASSP3P3S1C "assp3p3s1c"
810: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
811: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
812: #define TSROSWARK3 "ark3"
813: #define TSROSWTHETA1 "theta1"
814: #define TSROSWTHETA2 "theta2"
815: #define TSROSWGRK4T "grk4t"
816: #define TSROSWSHAMP4 "shamp4"
817: #define TSROSWVELDD4 "veldd4"
818: #define TSROSW4L "4l"
820: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
821: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
822: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
823: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
824: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
825: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
826: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
827: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
829: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS,PetscInt);
830: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS,PetscInt*);
832: /*
833: PETSc interface to Sundials
834: */
835: #ifdef PETSC_HAVE_SUNDIALS
836: typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
837: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
838: typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
839: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
840: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
841: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
842: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
843: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
844: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
845: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
846: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
847: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
848: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
849: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
850: PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
851: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
852: #endif
854: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
855: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
856: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
857: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
859: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
860: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
861: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
863: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS,PetscReal);
864: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS,PetscReal,PetscReal,PetscReal,PetscReal);
865: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
867: PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
868: PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
870: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
871: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
873: #endif