Actual source code: petscts.h
petsc-dev 2012-05-24
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
17: S*/
18: typedef struct _p_TS* TS;
20: /*J
21: TSType - String with the name of a PETSc TS method or the creation function
22: with an optional dynamic library name, for example
23: http://www.mcs.anl.gov/petsc/lib.a:mytscreate()
25: Level: beginner
27: .seealso: TSSetType(), TS
28: J*/
29: #define TSType char*
30: #define TSEULER "euler"
31: #define TSBEULER "beuler"
32: #define TSPSEUDO "pseudo"
33: #define TSCN "cn"
34: #define TSSUNDIALS "sundials"
35: #define TSRK "rk"
36: #define TSPYTHON "python"
37: #define TSTHETA "theta"
38: #define TSALPHA "alpha"
39: #define TSGL "gl"
40: #define TSSSP "ssp"
41: #define TSARKIMEX "arkimex"
42: #define TSROSW "rosw"
44: /*E
45: TSProblemType - Determines the type of problem this TS object is to be used to solve
47: Level: beginner
49: .seealso: TSCreate()
50: E*/
51: typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
53: /*E
54: TSConvergedReason - reason a TS method has converged or not
56: Level: beginner
58: Developer Notes: this must match finclude/petscts.h
60: Each reason has its own manual page.
62: .seealso: TSGetConvergedReason()
63: E*/
64: typedef enum {
65: TS_CONVERGED_ITERATING = 0,
66: TS_CONVERGED_TIME = 1,
67: TS_CONVERGED_ITS = 2,
68: TS_DIVERGED_NONLINEAR_SOLVE = -1,
69: TS_DIVERGED_STEP_REJECTED = -2
70: } TSConvergedReason;
71: PETSC_EXTERN const char *const*TSConvergedReasons;
73: /*MC
74: TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
76: Level: beginner
78: .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt()
79: M*/
81: /*MC
82: TS_CONVERGED_TIME - the final time was reached
84: Level: beginner
86: .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration()
87: M*/
89: /*MC
90: TS_CONVERGED_ITS - the maximum number of iterations was reached prior to the final time
92: Level: beginner
94: .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration()
95: M*/
97: /*MC
98: TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
100: Level: beginner
102: .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason()
103: M*/
105: /*MC
106: TS_DIVERGED_STEP_REJECTED - too many steps were rejected
108: Level: beginner
110: .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt()
111: M*/
113: /* Logging support */
114: PETSC_EXTERN PetscClassId TS_CLASSID;
116: PETSC_EXTERN PetscErrorCode TSInitializePackage(const char[]);
118: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
119: PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
121: PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
122: PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
123: PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
124: PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
125: PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
127: PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
128: PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
129: PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
130: PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
131: PETSC_EXTERN PetscErrorCode TSSetUp(TS);
132: PETSC_EXTERN PetscErrorCode TSReset(TS);
134: PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
135: PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
137: PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
138: PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
139: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,PetscBool);
141: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
142: PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
143: PETSC_EXTERN PetscErrorCode TSMonitorSolutionCreate(TS,PetscViewer,void**);
144: PETSC_EXTERN PetscErrorCode TSMonitorSolutionDestroy(void**);
145: PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
146: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
147: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
149: PETSC_EXTERN PetscErrorCode TSStep(TS);
150: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
151: PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec,PetscReal*);
152: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
153: PETSC_EXTERN PetscErrorCode TSGetNonlinearSolveIterations(TS,PetscInt*);
154: PETSC_EXTERN PetscErrorCode TSGetLinearSolveIterations(TS,PetscInt*);
156: PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
157: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
158: PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
159: PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
160: PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
161: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
163: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
164: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
165: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
166: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
167: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
168: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
170: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
171: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
172: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
173: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
174: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
175: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
177: PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
178: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
179: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
180: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
182: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
183: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
184: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
185: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
186: PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
187: PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
188: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
189: PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*);
190: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
191: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
193: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
194: PETSC_EXTERN PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
195: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
196: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
198: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
199: PETSC_EXTERN PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
200: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
201: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
202: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
204: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
206: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
207: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
208: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
209: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool);
211: PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
213: /* Dynamic creation and loading functions */
214: PETSC_EXTERN PetscFList TSList;
215: PETSC_EXTERN PetscBool TSRegisterAllCalled;
216: PETSC_EXTERN PetscErrorCode TSGetType(TS,const TSType*);
217: PETSC_EXTERN PetscErrorCode TSSetType(TS,const TSType);
218: PETSC_EXTERN PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
219: PETSC_EXTERN PetscErrorCode TSRegisterAll(const char[]);
220: PETSC_EXTERN PetscErrorCode TSRegisterDestroy(void);
222: /*MC
223: TSRegisterDynamic - Adds a creation method to the TS package.
225: Synopsis:
226: PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
228: Not Collective
230: Input Parameters:
231: + name - The name of a new user-defined creation routine
232: . path - The path (either absolute or relative) of the library containing this routine
233: . func_name - The name of the creation routine
234: - create_func - The creation routine itself
236: Notes:
237: TSRegisterDynamic() may be called multiple times to add several user-defined tses.
239: If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
241: Sample usage:
242: .vb
243: TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
244: .ve
246: Then, your ts type can be chosen with the procedural interface via
247: .vb
248: TS ts;
249: TSCreate(MPI_Comm, &ts);
250: TSSetType(ts, "my_ts")
251: .ve
252: or at runtime via the option
253: .vb
254: -ts_type my_ts
255: .ve
257: Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
258: If your function is not being put into a shared library then use TSRegister() instead
260: Level: advanced
262: .keywords: TS, register
263: .seealso: TSRegisterAll(), TSRegisterDestroy()
264: M*/
265: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
266: #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
267: #else
268: #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
269: #endif
271: PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
272: PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
274: PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
276: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
277: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
279: PETSC_EXTERN PetscErrorCode TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
280: PETSC_EXTERN PetscErrorCode TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
281: PETSC_EXTERN PetscErrorCode TSMonitorLGDestroy(PetscDrawLG*);
283: /*J
284: TSSSPType - string with the name of TSSSP scheme.
286: Level: beginner
288: .seealso: TSSSPSetType(), TS
289: J*/
290: #define TSSSPType char*
291: #define TSSSPRKS2 "rks2"
292: #define TSSSPRKS3 "rks3"
293: #define TSSSPRK104 "rk104"
295: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,const TSSSPType);
296: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,const TSSSPType*);
297: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
298: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
300: /*S
301: TSAdapt - Abstract object that manages time-step adaptivity
303: Level: beginner
305: .seealso: TS, TSAdaptCreate(), TSAdaptType
306: S*/
307: typedef struct _p_TSAdapt *TSAdapt;
309: /*E
310: TSAdaptType - String with the name of TSAdapt scheme or the creation function
311: with an optional dynamic library name, for example
312: http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
314: Level: beginner
316: .seealso: TSAdaptSetType(), TS
317: E*/
318: #define TSAdaptType char*
319: #define TSADAPTBASIC "basic"
320: #define TSADAPTNONE "none"
321: #define TSADAPTCFL "cfl"
323: /*MC
324: TSAdaptRegisterDynamic - adds a TSAdapt implementation
326: Synopsis:
327: PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
329: Not Collective
331: Input Parameters:
332: + name_scheme - name of user-defined adaptivity scheme
333: . path - path (either absolute or relative) the library containing this scheme
334: . name_create - name of routine to create method context
335: - routine_create - routine to create method context
337: Notes:
338: TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
340: If dynamic libraries are used, then the fourth input argument (routine_create)
341: is ignored.
343: Sample usage:
344: .vb
345: TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
346: "MySchemeCreate",MySchemeCreate);
347: .ve
349: Then, your scheme can be chosen with the procedural interface via
350: $ TSAdaptSetType(ts,"my_scheme")
351: or at runtime via the option
352: $ -ts_adapt_type my_scheme
354: Level: advanced
356: Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
357: and others of the form ${any_environmental_variable} occuring in pathname will be
358: replaced with appropriate values.
360: .keywords: TSAdapt, register
362: .seealso: TSAdaptRegisterAll()
363: M*/
364: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
365: # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,0)
366: #else
367: # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,d)
368: #endif
370: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
371: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt));
372: PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(const char[]);
373: PETSC_EXTERN PetscErrorCode TSAdaptRegisterDestroy(void);
374: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(const char[]);
375: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
376: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
377: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,const TSAdaptType);
378: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
379: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
380: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
381: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
382: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
383: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
384: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
385: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt);
386: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
387: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
388: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
390: /*S
391: TSGLAdapt - Abstract object that manages time-step adaptivity
393: Level: beginner
395: Developer Notes:
396: This functionality should be replaced by the TSAdapt.
398: .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
399: S*/
400: typedef struct _p_TSGLAdapt *TSGLAdapt;
402: /*J
403: TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
404: with an optional dynamic library name, for example
405: http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
407: Level: beginner
409: .seealso: TSGLAdaptSetType(), TS
410: J*/
411: #define TSGLAdaptType char*
412: #define TSGLADAPT_NONE "none"
413: #define TSGLADAPT_SIZE "size"
414: #define TSGLADAPT_BOTH "both"
416: /*MC
417: TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
419: Synopsis:
420: PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
422: Not Collective
424: Input Parameters:
425: + name_scheme - name of user-defined adaptivity scheme
426: . path - path (either absolute or relative) the library containing this scheme
427: . name_create - name of routine to create method context
428: - routine_create - routine to create method context
430: Notes:
431: TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
433: If dynamic libraries are used, then the fourth input argument (routine_create)
434: is ignored.
436: Sample usage:
437: .vb
438: TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
439: "MySchemeCreate",MySchemeCreate);
440: .ve
442: Then, your scheme can be chosen with the procedural interface via
443: $ TSGLAdaptSetType(ts,"my_scheme")
444: or at runtime via the option
445: $ -ts_adapt_type my_scheme
447: Level: advanced
449: Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
450: and others of the form ${any_environmental_variable} occuring in pathname will be
451: replaced with appropriate values.
453: .keywords: TSGLAdapt, register
455: .seealso: TSGLAdaptRegisterAll()
456: M*/
457: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
458: # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0)
459: #else
460: # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d)
461: #endif
463: PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
464: PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(const char[]);
465: PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterDestroy(void);
466: PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(const char[]);
467: PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
468: PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
469: PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
470: PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
471: PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
472: PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
473: PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt);
474: PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
476: /*J
477: TSGLAcceptType - String with the name of TSGLAccept scheme or the function
478: with an optional dynamic library name, for example
479: http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
481: Level: beginner
483: .seealso: TSGLSetAcceptType(), TS
484: J*/
485: #define TSGLAcceptType char*
486: #define TSGLACCEPT_ALWAYS "always"
488: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
489: PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
491: /*MC
492: TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
494: Synopsis:
495: PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
497: Not Collective
499: Input Parameters:
500: + name_scheme - name of user-defined acceptance scheme
501: . path - path (either absolute or relative) the library containing this scheme
502: . name_create - name of routine to create method context
503: - routine_create - routine to create method context
505: Notes:
506: TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
508: If dynamic libraries are used, then the fourth input argument (routine_create)
509: is ignored.
511: Sample usage:
512: .vb
513: TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
514: "MySchemeCreate",MySchemeCreate);
515: .ve
517: Then, your scheme can be chosen with the procedural interface via
518: $ TSGLSetAcceptType(ts,"my_scheme")
519: or at runtime via the option
520: $ -ts_gl_accept_type my_scheme
522: Level: advanced
524: Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
525: and others of the form ${any_environmental_variable} occuring in pathname will be
526: replaced with appropriate values.
528: .keywords: TSGL, TSGLAcceptType, register
530: .seealso: TSGLRegisterAll()
531: M*/
532: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
533: # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
534: #else
535: # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
536: #endif
538: /*J
539: TSGLType - family of time integration method within the General Linear class
541: Level: beginner
543: .seealso: TSGLSetType(), TSGLRegister()
544: J*/
545: #define TSGLType char*
546: #define TSGL_IRKS "irks"
548: /*MC
549: TSGLRegisterDynamic - adds a TSGL implementation
551: Synopsis:
552: PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
554: Not Collective
556: Input Parameters:
557: + name_scheme - name of user-defined general linear scheme
558: . path - path (either absolute or relative) the library containing this scheme
559: . name_create - name of routine to create method context
560: - routine_create - routine to create method context
562: Notes:
563: TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
565: If dynamic libraries are used, then the fourth input argument (routine_create)
566: is ignored.
568: Sample usage:
569: .vb
570: TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
571: "MySchemeCreate",MySchemeCreate);
572: .ve
574: Then, your scheme can be chosen with the procedural interface via
575: $ TSGLSetType(ts,"my_scheme")
576: or at runtime via the option
577: $ -ts_gl_type my_scheme
579: Level: advanced
581: Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
582: and others of the form ${any_environmental_variable} occuring in pathname will be
583: replaced with appropriate values.
585: .keywords: TSGL, register
587: .seealso: TSGLRegisterAll()
588: M*/
589: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
590: # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0)
591: #else
592: # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d)
593: #endif
595: PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
596: PETSC_EXTERN PetscErrorCode TSGLRegisterAll(const char[]);
597: PETSC_EXTERN PetscErrorCode TSGLRegisterDestroy(void);
598: PETSC_EXTERN PetscErrorCode TSGLInitializePackage(const char[]);
599: PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
600: PETSC_EXTERN PetscErrorCode TSGLSetType(TS,const TSGLType);
601: PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
602: PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,const TSGLAcceptType);
604: /*J
605: TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
607: Level: beginner
609: .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
610: J*/
611: #define TSARKIMEXType char*
612: #define TSARKIMEXA2 "a2"
613: #define TSARKIMEXL2 "l2"
614: #define TSARKIMEXARS122 "ars122"
615: #define TSARKIMEX2C "2c"
616: #define TSARKIMEX2D "2d"
617: #define TSARKIMEX2E "2e"
618: #define TSARKIMEXPRSSP2 "prssp2"
619: #define TSARKIMEX3 "3"
620: #define TSARKIMEXBPR3 "bpr3"
621: #define TSARKIMEXARS443 "ars443"
622: #define TSARKIMEX4 "4"
623: #define TSARKIMEX5 "5"
624: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
625: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
626: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
627: PETSC_EXTERN PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
628: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
629: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
630: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
631: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void);
633: /*J
634: TSRosWType - String with the name of a Rosenbrock-W method.
636: Level: beginner
638: .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
639: J*/
640: #define TSRosWType char*
641: #define TSROSW2M "2m"
642: #define TSROSW2P "2p"
643: #define TSROSWRA3PW "ra3pw"
644: #define TSROSWRA34PW2 "ra34pw2"
645: #define TSROSWRODAS3 "rodas3"
646: #define TSROSWSANDU3 "sandu3"
647: #define TSROSWASSP3P3S1C "assp3p3s1c"
648: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
649: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
650: #define TSROSWARK3 "ark3"
651: #define TSROSWTHETA1 "theta1"
652: #define TSROSWTHETA2 "theta2"
655: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,const TSRosWType*);
656: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,const TSRosWType);
657: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
658: PETSC_EXTERN PetscErrorCode TSRosWRegister(const TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
659: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
660: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(const char path[]);
661: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
662: PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
664: /*
665: PETSc interface to Sundials
666: */
667: #ifdef PETSC_HAVE_SUNDIALS
668: typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
669: PETSC_EXTERN const char *TSSundialsLmmTypes[];
670: typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
671: PETSC_EXTERN const char *TSSundialsGramSchmidtTypes[];
672: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
673: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
674: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
675: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
676: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
677: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
678: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
679: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
680: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
681: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
682: PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
683: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
684: #endif
686: PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal);
688: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
689: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
690: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
691: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
693: PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
694: PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
695: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
696: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
697: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
699: PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
700: PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
702: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
703: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
705: #endif