xref: /petsc/src/ts/interface/ts.c (revision 1e3347e8910e9540ca4e5e63802e011468fd198c)
163dd3a1aSKris Buschelman #define PETSCTS_DLL
263dd3a1aSKris Buschelman 
3e090d566SSatish Balay #include "src/ts/tsimpl.h"        /*I "petscts.h"  I*/
4d763cef2SBarry Smith 
5d5ba7fb7SMatthew Knepley /* Logging support */
663dd3a1aSKris Buschelman PetscCookie PETSCTS_DLLEXPORT TS_COOKIE = 0;
76849ba73SBarry Smith PetscEvent  TS_Step = 0, TS_PseudoComputeTimeStep = 0, TS_FunctionEval = 0, TS_JacobianEval = 0;
8d405a339SMatthew Knepley 
94a2ae208SSatish Balay #undef __FUNCT__
10bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
11bdad233fSMatthew Knepley /*
12bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
13bdad233fSMatthew Knepley 
14bdad233fSMatthew Knepley   Collective on TS
15bdad233fSMatthew Knepley 
16bdad233fSMatthew Knepley   Input Parameter:
17bdad233fSMatthew Knepley . ts - The ts
18bdad233fSMatthew Knepley 
19bdad233fSMatthew Knepley   Level: intermediate
20bdad233fSMatthew Knepley 
21bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
22bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
23bdad233fSMatthew Knepley */
246849ba73SBarry Smith static PetscErrorCode TSSetTypeFromOptions(TS ts)
25bdad233fSMatthew Knepley {
26bdad233fSMatthew Knepley   PetscTruth     opt;
272fc52814SBarry Smith   const char     *defaultType;
28bdad233fSMatthew Knepley   char           typeName[256];
29dfbe8321SBarry Smith   PetscErrorCode ierr;
30bdad233fSMatthew Knepley 
31bdad233fSMatthew Knepley   PetscFunctionBegin;
32abc0a331SBarry Smith   if (ts->type_name) {
33bdad233fSMatthew Knepley     defaultType = ts->type_name;
34bdad233fSMatthew Knepley   } else {
35bdad233fSMatthew Knepley     defaultType = TS_EULER;
36bdad233fSMatthew Knepley   }
37bdad233fSMatthew Knepley 
38bdad233fSMatthew Knepley   if (!TSRegisterAllCalled) {
39bdad233fSMatthew Knepley     ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);
40bdad233fSMatthew Knepley   }
41bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
42a7cc72afSBarry Smith   if (opt) {
43bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
44bdad233fSMatthew Knepley   } else {
45bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
46bdad233fSMatthew Knepley   }
47bdad233fSMatthew Knepley   PetscFunctionReturn(0);
48bdad233fSMatthew Knepley }
49bdad233fSMatthew Knepley 
50bdad233fSMatthew Knepley #undef __FUNCT__
51bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
52bdad233fSMatthew Knepley /*@
53bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
54bdad233fSMatthew Knepley 
55bdad233fSMatthew Knepley    Collective on TS
56bdad233fSMatthew Knepley 
57bdad233fSMatthew Knepley    Input Parameter:
58bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
59bdad233fSMatthew Knepley 
60bdad233fSMatthew Knepley    Options Database Keys:
610f3b3ca1SHong Zhang +  -ts_type <type> - TS_EULER, TS_BEULER, TS_SUNDIALS, TS_PSEUDO, TS_CRANK_NICHOLSON
62bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
63bdad233fSMatthew Knepley .  -ts_max_time time - maximum time to compute to
64bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
65bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
66bdad233fSMatthew Knepley -  -ts_xmonitor - plot information at each timestep
67bdad233fSMatthew Knepley 
68bdad233fSMatthew Knepley    Level: beginner
69bdad233fSMatthew Knepley 
70bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
71bdad233fSMatthew Knepley 
72bdad233fSMatthew Knepley .seealso: TSGetType
73bdad233fSMatthew Knepley @*/
7463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetFromOptions(TS ts)
75bdad233fSMatthew Knepley {
76bdad233fSMatthew Knepley   PetscReal      dt;
77eabae89aSBarry Smith   PetscTruth     opt,flg;
78dfbe8321SBarry Smith   PetscErrorCode ierr;
79eabae89aSBarry Smith   PetscViewer    monviewer;
80eabae89aSBarry Smith   char           monfilename[PETSC_MAX_PATH_LEN];
81bdad233fSMatthew Knepley 
82bdad233fSMatthew Knepley   PetscFunctionBegin;
834482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
84bdad233fSMatthew Knepley   ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");CHKERRQ(ierr);
85bdad233fSMatthew Knepley 
86bdad233fSMatthew Knepley     /* Handle generic TS options */
87bdad233fSMatthew Knepley     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
88bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
89bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
90bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
91a7cc72afSBarry Smith     if (opt) {
92bdad233fSMatthew Knepley       ts->initial_time_step = ts->time_step = dt;
93bdad233fSMatthew Knepley     }
94bdad233fSMatthew Knepley 
95bdad233fSMatthew Knepley     /* Monitor options */
96eabae89aSBarry Smith     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSDefaultMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
97eabae89aSBarry Smith     if (flg) {
98eabae89aSBarry Smith       ierr = PetscViewerASCIIOpen(ts->comm,monfilename,&monviewer);CHKERRQ(ierr);
99eabae89aSBarry Smith       ierr = TSSetMonitor(ts,TSDefaultMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
100bdad233fSMatthew Knepley     }
101bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_xmonitor","Monitor timestep size graphically","TSLGMonitor",&opt);CHKERRQ(ierr);
102a7cc72afSBarry Smith     if (opt) {
103bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
104bdad233fSMatthew Knepley     }
105bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_vecmonitor","Monitor solution graphically","TSVecViewMonitor",&opt);CHKERRQ(ierr);
106a7cc72afSBarry Smith     if (opt) {
107bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSVecViewMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
108bdad233fSMatthew Knepley     }
109bdad233fSMatthew Knepley 
110bdad233fSMatthew Knepley     /* Handle TS type options */
111bdad233fSMatthew Knepley     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
112bdad233fSMatthew Knepley 
113bdad233fSMatthew Knepley     /* Handle specific TS options */
114abc0a331SBarry Smith     if (ts->ops->setfromoptions) {
115bdad233fSMatthew Knepley       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
116bdad233fSMatthew Knepley     }
117bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
118bdad233fSMatthew Knepley 
119bdad233fSMatthew Knepley   /* Handle subobject options */
120bdad233fSMatthew Knepley   switch(ts->problem_type) {
121156fc9a6SMatthew Knepley     /* Should check for implicit/explicit */
122bdad233fSMatthew Knepley   case TS_LINEAR:
123abc0a331SBarry Smith     if (ts->ksp) {
1247c236d22SBarry Smith       ierr = KSPSetOperators(ts->ksp,ts->A,ts->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
12594b7f48cSBarry Smith       ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr);
126156fc9a6SMatthew Knepley     }
127bdad233fSMatthew Knepley     break;
128bdad233fSMatthew Knepley   case TS_NONLINEAR:
129abc0a331SBarry Smith     if (ts->snes) {
1307c236d22SBarry Smith       /* this is a bit of a hack, but it gets the matrix information into SNES earlier
1317c236d22SBarry Smith          so that SNES and KSP have more information to pick reasonable defaults
1327c236d22SBarry Smith          before they allow users to set options */
1337c236d22SBarry Smith       ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,0,ts);CHKERRQ(ierr);
134bdad233fSMatthew Knepley       ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
135156fc9a6SMatthew Knepley     }
136bdad233fSMatthew Knepley     break;
137bdad233fSMatthew Knepley   default:
13877431f27SBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type);
139bdad233fSMatthew Knepley   }
140bdad233fSMatthew Knepley 
141bdad233fSMatthew Knepley   PetscFunctionReturn(0);
142bdad233fSMatthew Knepley }
143bdad233fSMatthew Knepley 
144bdad233fSMatthew Knepley #undef  __FUNCT__
145bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
146bdad233fSMatthew Knepley /*@
147bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
148bdad233fSMatthew Knepley 
149bdad233fSMatthew Knepley   Collective on TS
150bdad233fSMatthew Knepley 
151bdad233fSMatthew Knepley   Input Parameter:
152bdad233fSMatthew Knepley . ts - The ts
153bdad233fSMatthew Knepley 
154bdad233fSMatthew Knepley   Level: intermediate
155bdad233fSMatthew Knepley 
156bdad233fSMatthew Knepley .keywords: TS, view, options, database
157bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
158bdad233fSMatthew Knepley @*/
15963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSViewFromOptions(TS ts,const char title[])
160bdad233fSMatthew Knepley {
161bdad233fSMatthew Knepley   PetscViewer    viewer;
162bdad233fSMatthew Knepley   PetscDraw      draw;
163bdad233fSMatthew Knepley   PetscTruth     opt;
164e10c49a3SBarry Smith   char           fileName[PETSC_MAX_PATH_LEN];
165dfbe8321SBarry Smith   PetscErrorCode ierr;
166bdad233fSMatthew Knepley 
167bdad233fSMatthew Knepley   PetscFunctionBegin;
168eabae89aSBarry Smith   ierr = PetscOptionsGetString(ts->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr);
169eabae89aSBarry Smith   if (opt && !PetscPreLoadingOn) {
170eabae89aSBarry Smith     ierr = PetscViewerASCIIOpen(ts->comm,fileName,&viewer);CHKERRQ(ierr);
171bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
172bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
173bdad233fSMatthew Knepley   }
174bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);CHKERRQ(ierr);
175a7cc72afSBarry Smith   if (opt) {
176bdad233fSMatthew Knepley     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
177bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
178a7cc72afSBarry Smith     if (title) {
1791836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
180bdad233fSMatthew Knepley     } else {
181bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject) ts);CHKERRQ(ierr);
1821836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, ts->name);CHKERRQ(ierr);
183bdad233fSMatthew Knepley     }
184bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
185bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
186bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
187bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
188bdad233fSMatthew Knepley   }
189bdad233fSMatthew Knepley   PetscFunctionReturn(0);
190bdad233fSMatthew Knepley }
191bdad233fSMatthew Knepley 
192bdad233fSMatthew Knepley #undef __FUNCT__
1934a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
194a7a1495cSBarry Smith /*@
1958c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
196a7a1495cSBarry Smith       set with TSSetRHSJacobian().
197a7a1495cSBarry Smith 
198a7a1495cSBarry Smith    Collective on TS and Vec
199a7a1495cSBarry Smith 
200a7a1495cSBarry Smith    Input Parameters:
201a7a1495cSBarry Smith +  ts - the SNES context
202a7a1495cSBarry Smith .  t - current timestep
203a7a1495cSBarry Smith -  x - input vector
204a7a1495cSBarry Smith 
205a7a1495cSBarry Smith    Output Parameters:
206a7a1495cSBarry Smith +  A - Jacobian matrix
207a7a1495cSBarry Smith .  B - optional preconditioning matrix
208a7a1495cSBarry Smith -  flag - flag indicating matrix structure
209a7a1495cSBarry Smith 
210a7a1495cSBarry Smith    Notes:
211a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
212a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
213a7a1495cSBarry Smith 
21494b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
215a7a1495cSBarry Smith    flag parameter.
216a7a1495cSBarry Smith 
217a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
218a7a1495cSBarry Smith 
219a7a1495cSBarry Smith    Level: developer
220a7a1495cSBarry Smith 
221a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
222a7a1495cSBarry Smith 
22394b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
224a7a1495cSBarry Smith @*/
22563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
226a7a1495cSBarry Smith {
227dfbe8321SBarry Smith   PetscErrorCode ierr;
228a7a1495cSBarry Smith 
229a7a1495cSBarry Smith   PetscFunctionBegin;
2304482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2314482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
232c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
233a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
23429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
235a7a1495cSBarry Smith   }
236000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
237d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
238a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
239a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
240000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
241a7a1495cSBarry Smith     PetscStackPop;
242d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
243a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2444482741eSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE,4);
2454482741eSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE,5);
246ef66eb69SBarry Smith   } else {
247ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249ef66eb69SBarry Smith     if (*A != *B) {
250ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252ef66eb69SBarry Smith     }
253ef66eb69SBarry Smith   }
254a7a1495cSBarry Smith   PetscFunctionReturn(0);
255a7a1495cSBarry Smith }
256a7a1495cSBarry Smith 
2574a2ae208SSatish Balay #undef __FUNCT__
2584a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
259d763cef2SBarry Smith /*
260d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
261d763cef2SBarry Smith 
262d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
263d763cef2SBarry Smith    this routine applies the matrix.
264d763cef2SBarry Smith */
265dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
266d763cef2SBarry Smith {
267dfbe8321SBarry Smith   PetscErrorCode ierr;
268d763cef2SBarry Smith 
269d763cef2SBarry Smith   PetscFunctionBegin;
2704482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2714482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
2724482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
273d763cef2SBarry Smith 
274d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
275000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
276d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
277000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
278d763cef2SBarry Smith     PetscStackPop;
2797c922b88SBarry Smith   } else {
280000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
281d763cef2SBarry Smith       MatStructure flg;
282d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
283000e7ae3SMatthew Knepley       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
284d763cef2SBarry Smith       PetscStackPop;
285d763cef2SBarry Smith     }
286d763cef2SBarry Smith     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
2877c922b88SBarry Smith   }
288d763cef2SBarry Smith 
289d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
290d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
291d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
292d763cef2SBarry Smith 
293d763cef2SBarry Smith   PetscFunctionReturn(0);
294d763cef2SBarry Smith }
295d763cef2SBarry Smith 
2964a2ae208SSatish Balay #undef __FUNCT__
2974a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
298d763cef2SBarry Smith /*@C
299d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
300d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
301d763cef2SBarry Smith 
302d763cef2SBarry Smith     Collective on TS
303d763cef2SBarry Smith 
304d763cef2SBarry Smith     Input Parameters:
305d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
306d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
307d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
308d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
309d763cef2SBarry Smith 
310d763cef2SBarry Smith     Calling sequence of func:
31187828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
312d763cef2SBarry Smith 
313d763cef2SBarry Smith +   t - current timestep
314d763cef2SBarry Smith .   u - input vector
315d763cef2SBarry Smith .   F - function vector
316d763cef2SBarry Smith -   ctx - [optional] user-defined function context
317d763cef2SBarry Smith 
318d763cef2SBarry Smith     Important:
319d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
320d763cef2SBarry Smith 
321d763cef2SBarry Smith     Level: beginner
322d763cef2SBarry Smith 
323d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
324d763cef2SBarry Smith 
325d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
326d763cef2SBarry Smith @*/
32763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
328d763cef2SBarry Smith {
329d763cef2SBarry Smith   PetscFunctionBegin;
330d763cef2SBarry Smith 
3314482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
332d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
33329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
334d763cef2SBarry Smith   }
335000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
336d763cef2SBarry Smith   ts->funP             = ctx;
337d763cef2SBarry Smith   PetscFunctionReturn(0);
338d763cef2SBarry Smith }
339d763cef2SBarry Smith 
3404a2ae208SSatish Balay #undef __FUNCT__
3414a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
342d763cef2SBarry Smith /*@C
343d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
344d763cef2SBarry Smith    Also sets the location to store A.
345d763cef2SBarry Smith 
346d763cef2SBarry Smith    Collective on TS
347d763cef2SBarry Smith 
348d763cef2SBarry Smith    Input Parameters:
349d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
350d763cef2SBarry Smith .  A   - matrix
351d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
352d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
353d763cef2SBarry Smith          if A is not a function of t.
354d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
355d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
356d763cef2SBarry Smith 
357d763cef2SBarry Smith    Calling sequence of func:
358a7cc72afSBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
359d763cef2SBarry Smith 
360d763cef2SBarry Smith +  t - current timestep
361d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
362d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
363d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
36494b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
365d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
366d763cef2SBarry Smith 
367d763cef2SBarry Smith    Notes:
36894b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
369d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
370d763cef2SBarry Smith 
371d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
372d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
373d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
374d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
375d763cef2SBarry Smith    throughout the global iterations.
376d763cef2SBarry Smith 
377d763cef2SBarry Smith    Important:
378d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
379d763cef2SBarry Smith 
380d763cef2SBarry Smith    Level: beginner
381d763cef2SBarry Smith 
382d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
383d763cef2SBarry Smith 
384d763cef2SBarry Smith .seealso: TSSetRHSFunction()
385d763cef2SBarry Smith @*/
38663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSMatrix(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
387d763cef2SBarry Smith {
388d763cef2SBarry Smith   PetscFunctionBegin;
3894482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
3904482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
3914482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
392c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
393c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
394d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
39529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
396d763cef2SBarry Smith   }
397d763cef2SBarry Smith 
398000e7ae3SMatthew Knepley   ts->ops->rhsmatrix = f;
399d763cef2SBarry Smith   ts->jacP           = ctx;
400d763cef2SBarry Smith   ts->A              = A;
401d763cef2SBarry Smith   ts->B              = B;
402d763cef2SBarry Smith 
403d763cef2SBarry Smith   PetscFunctionReturn(0);
404d763cef2SBarry Smith }
405d763cef2SBarry Smith 
4064a2ae208SSatish Balay #undef __FUNCT__
4074a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
408d763cef2SBarry Smith /*@C
409d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
410d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
4113c94ec11SBarry Smith    Use TSSetRHSMatrix() for linear problems.
412d763cef2SBarry Smith 
413d763cef2SBarry Smith    Collective on TS
414d763cef2SBarry Smith 
415d763cef2SBarry Smith    Input Parameters:
416d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
417d763cef2SBarry Smith .  A   - Jacobian matrix
418d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
419d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
420d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
421d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
422d763cef2SBarry Smith 
423d763cef2SBarry Smith    Calling sequence of func:
42487828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
425d763cef2SBarry Smith 
426d763cef2SBarry Smith +  t - current timestep
427d763cef2SBarry Smith .  u - input vector
428d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
429d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
430d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
43194b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
432d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
433d763cef2SBarry Smith 
434d763cef2SBarry Smith    Notes:
43594b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
436d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
437d763cef2SBarry Smith 
438d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
439d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
440d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
441d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
442d763cef2SBarry Smith    throughout the global iterations.
443d763cef2SBarry Smith 
444d763cef2SBarry Smith    Level: beginner
445d763cef2SBarry Smith 
446d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
447d763cef2SBarry Smith 
448d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
4493c94ec11SBarry Smith           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetRHSMatrix()
450d763cef2SBarry Smith 
451d763cef2SBarry Smith @*/
45263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
453d763cef2SBarry Smith {
454d763cef2SBarry Smith   PetscFunctionBegin;
4554482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4564482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
4574482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
458c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
459c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
460d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
46129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
462d763cef2SBarry Smith   }
463d763cef2SBarry Smith 
464000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
465d763cef2SBarry Smith   ts->jacP             = ctx;
466d763cef2SBarry Smith   ts->A                = A;
467d763cef2SBarry Smith   ts->B                = B;
468d763cef2SBarry Smith   PetscFunctionReturn(0);
469d763cef2SBarry Smith }
470d763cef2SBarry Smith 
4714a2ae208SSatish Balay #undef __FUNCT__
4724a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
473d763cef2SBarry Smith /*
474d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
475d763cef2SBarry Smith 
476d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
477d763cef2SBarry Smith    this routine applies the matrix.
478d763cef2SBarry Smith */
479dfbe8321SBarry Smith PetscErrorCode TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
480d763cef2SBarry Smith {
481dfbe8321SBarry Smith   PetscErrorCode ierr;
482d763cef2SBarry Smith 
483d763cef2SBarry Smith   PetscFunctionBegin;
4844482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4854482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
486c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,x,3);
487d763cef2SBarry Smith 
488000e7ae3SMatthew Knepley   if (ts->ops->rhsbc) {
489d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
490000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
491d763cef2SBarry Smith     PetscStackPop;
492d763cef2SBarry Smith     PetscFunctionReturn(0);
493d763cef2SBarry Smith   }
494d763cef2SBarry Smith 
495d763cef2SBarry Smith   PetscFunctionReturn(0);
496d763cef2SBarry Smith }
497d763cef2SBarry Smith 
4984a2ae208SSatish Balay #undef __FUNCT__
4994a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
500d763cef2SBarry Smith /*@C
501d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
502d763cef2SBarry Smith     boundary conditions for the function F.
503d763cef2SBarry Smith 
504d763cef2SBarry Smith     Collective on TS
505d763cef2SBarry Smith 
506d763cef2SBarry Smith     Input Parameters:
507d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
508d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
509d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
510d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
511d763cef2SBarry Smith 
512d763cef2SBarry Smith     Calling sequence of func:
51387828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
514d763cef2SBarry Smith 
515d763cef2SBarry Smith +   t - current timestep
516d763cef2SBarry Smith .   F - function vector
517d763cef2SBarry Smith -   ctx - [optional] user-defined function context
518d763cef2SBarry Smith 
519d763cef2SBarry Smith     Level: intermediate
520d763cef2SBarry Smith 
521d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
522d763cef2SBarry Smith @*/
52363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSBoundaryConditions(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
524d763cef2SBarry Smith {
525d763cef2SBarry Smith   PetscFunctionBegin;
526d763cef2SBarry Smith 
5274482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
528d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
52929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
530d763cef2SBarry Smith   }
531000e7ae3SMatthew Knepley   ts->ops->rhsbc = f;
532d763cef2SBarry Smith   ts->bcP        = ctx;
533d763cef2SBarry Smith   PetscFunctionReturn(0);
534d763cef2SBarry Smith }
535d763cef2SBarry Smith 
5364a2ae208SSatish Balay #undef __FUNCT__
5374a2ae208SSatish Balay #define __FUNCT__ "TSView"
5387e2c5f70SBarry Smith /*@C
539d763cef2SBarry Smith     TSView - Prints the TS data structure.
540d763cef2SBarry Smith 
5414c49b128SBarry Smith     Collective on TS
542d763cef2SBarry Smith 
543d763cef2SBarry Smith     Input Parameters:
544d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
545d763cef2SBarry Smith -   viewer - visualization context
546d763cef2SBarry Smith 
547d763cef2SBarry Smith     Options Database Key:
548d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
549d763cef2SBarry Smith 
550d763cef2SBarry Smith     Notes:
551d763cef2SBarry Smith     The available visualization contexts include
552b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
553b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
554d763cef2SBarry Smith          output where only the first processor opens
555d763cef2SBarry Smith          the file.  All other processors send their
556d763cef2SBarry Smith          data to the first processor to print.
557d763cef2SBarry Smith 
558d763cef2SBarry Smith     The user can open an alternative visualization context with
559b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
560d763cef2SBarry Smith 
561d763cef2SBarry Smith     Level: beginner
562d763cef2SBarry Smith 
563d763cef2SBarry Smith .keywords: TS, timestep, view
564d763cef2SBarry Smith 
565b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
566d763cef2SBarry Smith @*/
56763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer)
568d763cef2SBarry Smith {
569dfbe8321SBarry Smith   PetscErrorCode ierr;
570454a90a3SBarry Smith   char           *type;
57132077d6dSBarry Smith   PetscTruth     iascii,isstring;
572d763cef2SBarry Smith 
573d763cef2SBarry Smith   PetscFunctionBegin;
5744482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
575b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
5764482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
577c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
578fd16b177SBarry Smith 
57932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
580b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
58132077d6dSBarry Smith   if (iascii) {
582b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
583454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
584454a90a3SBarry Smith     if (type) {
585b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
586184914b5SBarry Smith     } else {
587b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
588184914b5SBarry Smith     }
589000e7ae3SMatthew Knepley     if (ts->ops->view) {
590b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
591000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
592b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
593d763cef2SBarry Smith     }
59477431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
595a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
596d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
59777431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
598d763cef2SBarry Smith     }
59977431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
6000f5bd95cSBarry Smith   } else if (isstring) {
601454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
602b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
603d763cef2SBarry Smith   }
604b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
60594b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
606d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
607b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
608d763cef2SBarry Smith   PetscFunctionReturn(0);
609d763cef2SBarry Smith }
610d763cef2SBarry Smith 
611d763cef2SBarry Smith 
6124a2ae208SSatish Balay #undef __FUNCT__
6134a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
614d763cef2SBarry Smith /*@C
615d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
616d763cef2SBarry Smith    the timesteppers.
617d763cef2SBarry Smith 
618d763cef2SBarry Smith    Collective on TS
619d763cef2SBarry Smith 
620d763cef2SBarry Smith    Input Parameters:
621d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
622d763cef2SBarry Smith -  usrP - optional user context
623d763cef2SBarry Smith 
624d763cef2SBarry Smith    Level: intermediate
625d763cef2SBarry Smith 
626d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
627d763cef2SBarry Smith 
628d763cef2SBarry Smith .seealso: TSGetApplicationContext()
629d763cef2SBarry Smith @*/
63063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP)
631d763cef2SBarry Smith {
632d763cef2SBarry Smith   PetscFunctionBegin;
6334482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
634d763cef2SBarry Smith   ts->user = usrP;
635d763cef2SBarry Smith   PetscFunctionReturn(0);
636d763cef2SBarry Smith }
637d763cef2SBarry Smith 
6384a2ae208SSatish Balay #undef __FUNCT__
6394a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
640d763cef2SBarry Smith /*@C
641d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
642d763cef2SBarry Smith     timestepper.
643d763cef2SBarry Smith 
644d763cef2SBarry Smith     Not Collective
645d763cef2SBarry Smith 
646d763cef2SBarry Smith     Input Parameter:
647d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
648d763cef2SBarry Smith 
649d763cef2SBarry Smith     Output Parameter:
650d763cef2SBarry Smith .   usrP - user context
651d763cef2SBarry Smith 
652d763cef2SBarry Smith     Level: intermediate
653d763cef2SBarry Smith 
654d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
655d763cef2SBarry Smith 
656d763cef2SBarry Smith .seealso: TSSetApplicationContext()
657d763cef2SBarry Smith @*/
65863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP)
659d763cef2SBarry Smith {
660d763cef2SBarry Smith   PetscFunctionBegin;
6614482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
662d763cef2SBarry Smith   *usrP = ts->user;
663d763cef2SBarry Smith   PetscFunctionReturn(0);
664d763cef2SBarry Smith }
665d763cef2SBarry Smith 
6664a2ae208SSatish Balay #undef __FUNCT__
6674a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
668d763cef2SBarry Smith /*@
669d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
670d763cef2SBarry Smith 
671d763cef2SBarry Smith    Not Collective
672d763cef2SBarry Smith 
673d763cef2SBarry Smith    Input Parameter:
674d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
675d763cef2SBarry Smith 
676d763cef2SBarry Smith    Output Parameter:
677d763cef2SBarry Smith .  iter - number steps so far
678d763cef2SBarry Smith 
679d763cef2SBarry Smith    Level: intermediate
680d763cef2SBarry Smith 
681d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
682d763cef2SBarry Smith @*/
68363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter)
684d763cef2SBarry Smith {
685d763cef2SBarry Smith   PetscFunctionBegin;
6864482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
6874482741eSBarry Smith   PetscValidIntPointer(iter,2);
688d763cef2SBarry Smith   *iter = ts->steps;
689d763cef2SBarry Smith   PetscFunctionReturn(0);
690d763cef2SBarry Smith }
691d763cef2SBarry Smith 
6924a2ae208SSatish Balay #undef __FUNCT__
6934a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
694d763cef2SBarry Smith /*@
695d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
696d763cef2SBarry Smith    as well as the initial time.
697d763cef2SBarry Smith 
698d763cef2SBarry Smith    Collective on TS
699d763cef2SBarry Smith 
700d763cef2SBarry Smith    Input Parameters:
701d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
702d763cef2SBarry Smith .  initial_time - the initial time
703d763cef2SBarry Smith -  time_step - the size of the timestep
704d763cef2SBarry Smith 
705d763cef2SBarry Smith    Level: intermediate
706d763cef2SBarry Smith 
707d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
708d763cef2SBarry Smith 
709d763cef2SBarry Smith .keywords: TS, set, initial, timestep
710d763cef2SBarry Smith @*/
71163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
712d763cef2SBarry Smith {
713d763cef2SBarry Smith   PetscFunctionBegin;
7144482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
715d763cef2SBarry Smith   ts->time_step         = time_step;
716d763cef2SBarry Smith   ts->initial_time_step = time_step;
717d763cef2SBarry Smith   ts->ptime             = initial_time;
718d763cef2SBarry Smith   PetscFunctionReturn(0);
719d763cef2SBarry Smith }
720d763cef2SBarry Smith 
7214a2ae208SSatish Balay #undef __FUNCT__
7224a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
723d763cef2SBarry Smith /*@
724d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
725d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
726d763cef2SBarry Smith 
727d763cef2SBarry Smith    Collective on TS
728d763cef2SBarry Smith 
729d763cef2SBarry Smith    Input Parameters:
730d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
731d763cef2SBarry Smith -  time_step - the size of the timestep
732d763cef2SBarry Smith 
733d763cef2SBarry Smith    Level: intermediate
734d763cef2SBarry Smith 
735d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
736d763cef2SBarry Smith 
737d763cef2SBarry Smith .keywords: TS, set, timestep
738d763cef2SBarry Smith @*/
73963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step)
740d763cef2SBarry Smith {
741d763cef2SBarry Smith   PetscFunctionBegin;
7424482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
743d763cef2SBarry Smith   ts->time_step = time_step;
744d763cef2SBarry Smith   PetscFunctionReturn(0);
745d763cef2SBarry Smith }
746d763cef2SBarry Smith 
7474a2ae208SSatish Balay #undef __FUNCT__
7484a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
749d763cef2SBarry Smith /*@
750d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
751d763cef2SBarry Smith 
752d763cef2SBarry Smith    Not Collective
753d763cef2SBarry Smith 
754d763cef2SBarry Smith    Input Parameter:
755d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
756d763cef2SBarry Smith 
757d763cef2SBarry Smith    Output Parameter:
758d763cef2SBarry Smith .  dt - the current timestep size
759d763cef2SBarry Smith 
760d763cef2SBarry Smith    Level: intermediate
761d763cef2SBarry Smith 
762d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
763d763cef2SBarry Smith 
764d763cef2SBarry Smith .keywords: TS, get, timestep
765d763cef2SBarry Smith @*/
76663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt)
767d763cef2SBarry Smith {
768d763cef2SBarry Smith   PetscFunctionBegin;
7694482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
7704482741eSBarry Smith   PetscValidDoublePointer(dt,2);
771d763cef2SBarry Smith   *dt = ts->time_step;
772d763cef2SBarry Smith   PetscFunctionReturn(0);
773d763cef2SBarry Smith }
774d763cef2SBarry Smith 
7754a2ae208SSatish Balay #undef __FUNCT__
7764a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
777d8e5e3e6SSatish Balay /*@
778d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
779d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
780d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
781d763cef2SBarry Smith    the solution at the next timestep has been calculated.
782d763cef2SBarry Smith 
783d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
784d763cef2SBarry Smith 
785d763cef2SBarry Smith    Input Parameter:
786d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
787d763cef2SBarry Smith 
788d763cef2SBarry Smith    Output Parameter:
789d763cef2SBarry Smith .  v - the vector containing the solution
790d763cef2SBarry Smith 
791d763cef2SBarry Smith    Level: intermediate
792d763cef2SBarry Smith 
793d763cef2SBarry Smith .seealso: TSGetTimeStep()
794d763cef2SBarry Smith 
795d763cef2SBarry Smith .keywords: TS, timestep, get, solution
796d763cef2SBarry Smith @*/
79763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v)
798d763cef2SBarry Smith {
799d763cef2SBarry Smith   PetscFunctionBegin;
8004482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
8014482741eSBarry Smith   PetscValidPointer(v,2);
802d763cef2SBarry Smith   *v = ts->vec_sol_always;
803d763cef2SBarry Smith   PetscFunctionReturn(0);
804d763cef2SBarry Smith }
805d763cef2SBarry Smith 
806bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
8074a2ae208SSatish Balay #undef __FUNCT__
808bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
809d8e5e3e6SSatish Balay /*@
810bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
811d763cef2SBarry Smith 
812bdad233fSMatthew Knepley   Not collective
813d763cef2SBarry Smith 
814bdad233fSMatthew Knepley   Input Parameters:
815bdad233fSMatthew Knepley + ts   - The TS
816bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
817d763cef2SBarry Smith .vb
818d763cef2SBarry Smith          U_t = A U
819d763cef2SBarry Smith          U_t = A(t) U
820d763cef2SBarry Smith          U_t = F(t,U)
821d763cef2SBarry Smith .ve
822d763cef2SBarry Smith 
823d763cef2SBarry Smith    Level: beginner
824d763cef2SBarry Smith 
825bdad233fSMatthew Knepley .keywords: TS, problem type
826bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
827d763cef2SBarry Smith @*/
82863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type)
829a7cc72afSBarry Smith {
830d763cef2SBarry Smith   PetscFunctionBegin;
8314482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
832bdad233fSMatthew Knepley   ts->problem_type = type;
833d763cef2SBarry Smith   PetscFunctionReturn(0);
834d763cef2SBarry Smith }
835d763cef2SBarry Smith 
836bdad233fSMatthew Knepley #undef __FUNCT__
837bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
838bdad233fSMatthew Knepley /*@C
839bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
840bdad233fSMatthew Knepley 
841bdad233fSMatthew Knepley   Not collective
842bdad233fSMatthew Knepley 
843bdad233fSMatthew Knepley   Input Parameter:
844bdad233fSMatthew Knepley . ts   - The TS
845bdad233fSMatthew Knepley 
846bdad233fSMatthew Knepley   Output Parameter:
847bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
848bdad233fSMatthew Knepley .vb
849bdad233fSMatthew Knepley          U_t = A U
850bdad233fSMatthew Knepley          U_t = A(t) U
851bdad233fSMatthew Knepley          U_t = F(t,U)
852bdad233fSMatthew Knepley .ve
853bdad233fSMatthew Knepley 
854bdad233fSMatthew Knepley    Level: beginner
855bdad233fSMatthew Knepley 
856bdad233fSMatthew Knepley .keywords: TS, problem type
857bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
858bdad233fSMatthew Knepley @*/
85963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type)
860a7cc72afSBarry Smith {
861bdad233fSMatthew Knepley   PetscFunctionBegin;
8624482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
8634482741eSBarry Smith   PetscValidIntPointer(type,2);
864bdad233fSMatthew Knepley   *type = ts->problem_type;
865bdad233fSMatthew Knepley   PetscFunctionReturn(0);
866bdad233fSMatthew Knepley }
867d763cef2SBarry Smith 
8684a2ae208SSatish Balay #undef __FUNCT__
8694a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
870d763cef2SBarry Smith /*@
871d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
872d763cef2SBarry Smith    of a timestepper.
873d763cef2SBarry Smith 
874d763cef2SBarry Smith    Collective on TS
875d763cef2SBarry Smith 
876d763cef2SBarry Smith    Input Parameter:
877d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
878d763cef2SBarry Smith 
879d763cef2SBarry Smith    Notes:
880d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
881d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
882d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
883d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
884d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
885d763cef2SBarry Smith 
886d763cef2SBarry Smith    Level: advanced
887d763cef2SBarry Smith 
888d763cef2SBarry Smith .keywords: TS, timestep, setup
889d763cef2SBarry Smith 
890d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
891d763cef2SBarry Smith @*/
89263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts)
893d763cef2SBarry Smith {
894dfbe8321SBarry Smith   PetscErrorCode ierr;
895d763cef2SBarry Smith 
896d763cef2SBarry Smith   PetscFunctionBegin;
8974482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
89829bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
899d763cef2SBarry Smith   if (!ts->type_name) {
900d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
901d763cef2SBarry Smith   }
902000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
903d763cef2SBarry Smith   ts->setupcalled = 1;
904d763cef2SBarry Smith   PetscFunctionReturn(0);
905d763cef2SBarry Smith }
906d763cef2SBarry Smith 
9074a2ae208SSatish Balay #undef __FUNCT__
9084a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
909d8e5e3e6SSatish Balay /*@
910d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
911d763cef2SBarry Smith    with TSCreate().
912d763cef2SBarry Smith 
913d763cef2SBarry Smith    Collective on TS
914d763cef2SBarry Smith 
915d763cef2SBarry Smith    Input Parameter:
916d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
917d763cef2SBarry Smith 
918d763cef2SBarry Smith    Level: beginner
919d763cef2SBarry Smith 
920d763cef2SBarry Smith .keywords: TS, timestepper, destroy
921d763cef2SBarry Smith 
922d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
923d763cef2SBarry Smith @*/
92463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts)
925d763cef2SBarry Smith {
9266849ba73SBarry Smith   PetscErrorCode ierr;
927d763cef2SBarry Smith 
928d763cef2SBarry Smith   PetscFunctionBegin;
9294482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
930d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
931d763cef2SBarry Smith 
932be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
9330f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
934be0abb6dSBarry Smith 
93594b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
936d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
937*1e3347e8SBarry Smith   if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);}
938d952e501SBarry Smith   ierr = TSClearMonitor(ts);CHKERRQ(ierr);
939a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
940d763cef2SBarry Smith   PetscFunctionReturn(0);
941d763cef2SBarry Smith }
942d763cef2SBarry Smith 
9434a2ae208SSatish Balay #undef __FUNCT__
9444a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
945d8e5e3e6SSatish Balay /*@
946d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
947d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
948d763cef2SBarry Smith 
949d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
950d763cef2SBarry Smith 
951d763cef2SBarry Smith    Input Parameter:
952d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
953d763cef2SBarry Smith 
954d763cef2SBarry Smith    Output Parameter:
955d763cef2SBarry Smith .  snes - the nonlinear solver context
956d763cef2SBarry Smith 
957d763cef2SBarry Smith    Notes:
958d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
959d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
96094b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
961d763cef2SBarry Smith 
962d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
963d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
964d763cef2SBarry Smith 
965d763cef2SBarry Smith    Level: beginner
966d763cef2SBarry Smith 
967d763cef2SBarry Smith .keywords: timestep, get, SNES
968d763cef2SBarry Smith @*/
96963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes)
970d763cef2SBarry Smith {
971d763cef2SBarry Smith   PetscFunctionBegin;
9724482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
9734482741eSBarry Smith   PetscValidPointer(snes,2);
97494b7f48cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
975d763cef2SBarry Smith   *snes = ts->snes;
976d763cef2SBarry Smith   PetscFunctionReturn(0);
977d763cef2SBarry Smith }
978d763cef2SBarry Smith 
9794a2ae208SSatish Balay #undef __FUNCT__
98094b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
981d8e5e3e6SSatish Balay /*@
98294b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
983d763cef2SBarry Smith    a TS (timestepper) context.
984d763cef2SBarry Smith 
98594b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
986d763cef2SBarry Smith 
987d763cef2SBarry Smith    Input Parameter:
988d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
989d763cef2SBarry Smith 
990d763cef2SBarry Smith    Output Parameter:
99194b7f48cSBarry Smith .  ksp - the nonlinear solver context
992d763cef2SBarry Smith 
993d763cef2SBarry Smith    Notes:
99494b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
995d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
996d763cef2SBarry Smith    KSP and PC contexts as well.
997d763cef2SBarry Smith 
99894b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
99994b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1000d763cef2SBarry Smith 
1001d763cef2SBarry Smith    Level: beginner
1002d763cef2SBarry Smith 
100394b7f48cSBarry Smith .keywords: timestep, get, KSP
1004d763cef2SBarry Smith @*/
100563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp)
1006d763cef2SBarry Smith {
1007d763cef2SBarry Smith   PetscFunctionBegin;
10084482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
10094482741eSBarry Smith   PetscValidPointer(ksp,2);
101029bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
101194b7f48cSBarry Smith   *ksp = ts->ksp;
1012d763cef2SBarry Smith   PetscFunctionReturn(0);
1013d763cef2SBarry Smith }
1014d763cef2SBarry Smith 
1015d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1016d763cef2SBarry Smith 
10174a2ae208SSatish Balay #undef __FUNCT__
1018adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1019adb62b0dSMatthew Knepley /*@
1020adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1021adb62b0dSMatthew Knepley    maximum time for iteration.
1022adb62b0dSMatthew Knepley 
1023adb62b0dSMatthew Knepley    Collective on TS
1024adb62b0dSMatthew Knepley 
1025adb62b0dSMatthew Knepley    Input Parameters:
1026adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1027adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1028adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1029adb62b0dSMatthew Knepley 
1030adb62b0dSMatthew Knepley    Level: intermediate
1031adb62b0dSMatthew Knepley 
1032adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1033adb62b0dSMatthew Knepley @*/
103463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1035adb62b0dSMatthew Knepley {
1036adb62b0dSMatthew Knepley   PetscFunctionBegin;
10374482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1038abc0a331SBarry Smith   if (maxsteps) {
10394482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1040adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1041adb62b0dSMatthew Knepley   }
1042abc0a331SBarry Smith   if (maxtime ) {
10434482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1044adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1045adb62b0dSMatthew Knepley   }
1046adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1047adb62b0dSMatthew Knepley }
1048adb62b0dSMatthew Knepley 
1049adb62b0dSMatthew Knepley #undef __FUNCT__
10504a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1051d763cef2SBarry Smith /*@
1052d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1053d763cef2SBarry Smith    maximum time for iteration.
1054d763cef2SBarry Smith 
1055d763cef2SBarry Smith    Collective on TS
1056d763cef2SBarry Smith 
1057d763cef2SBarry Smith    Input Parameters:
1058d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1059d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1060d763cef2SBarry Smith -  maxtime - final time to iterate to
1061d763cef2SBarry Smith 
1062d763cef2SBarry Smith    Options Database Keys:
1063d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1064d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1065d763cef2SBarry Smith 
1066d763cef2SBarry Smith    Notes:
1067d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1068d763cef2SBarry Smith 
1069d763cef2SBarry Smith    Level: intermediate
1070d763cef2SBarry Smith 
1071d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1072d763cef2SBarry Smith @*/
107363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1074d763cef2SBarry Smith {
1075d763cef2SBarry Smith   PetscFunctionBegin;
10764482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1077d763cef2SBarry Smith   ts->max_steps = maxsteps;
1078d763cef2SBarry Smith   ts->max_time  = maxtime;
1079d763cef2SBarry Smith   PetscFunctionReturn(0);
1080d763cef2SBarry Smith }
1081d763cef2SBarry Smith 
10824a2ae208SSatish Balay #undef __FUNCT__
10834a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1084d763cef2SBarry Smith /*@
1085d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1086d763cef2SBarry Smith    for use by the TS routines.
1087d763cef2SBarry Smith 
1088d763cef2SBarry Smith    Collective on TS and Vec
1089d763cef2SBarry Smith 
1090d763cef2SBarry Smith    Input Parameters:
1091d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1092d763cef2SBarry Smith -  x - the solution vector
1093d763cef2SBarry Smith 
1094d763cef2SBarry Smith    Level: beginner
1095d763cef2SBarry Smith 
1096d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1097d763cef2SBarry Smith @*/
109863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x)
1099d763cef2SBarry Smith {
1100d763cef2SBarry Smith   PetscFunctionBegin;
11014482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
11024482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1103d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1104d763cef2SBarry Smith   PetscFunctionReturn(0);
1105d763cef2SBarry Smith }
1106d763cef2SBarry Smith 
1107e74ef692SMatthew Knepley #undef __FUNCT__
1108e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultRhsBC"
1109000e7ae3SMatthew Knepley /*@
1110000e7ae3SMatthew Knepley   TSDefaultRhsBC - The default boundary condition function which does nothing.
1111000e7ae3SMatthew Knepley 
1112000e7ae3SMatthew Knepley   Collective on TS
1113000e7ae3SMatthew Knepley 
1114000e7ae3SMatthew Knepley   Input Parameters:
1115000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1116000e7ae3SMatthew Knepley . rhs - The Rhs
1117000e7ae3SMatthew Knepley - ctx - The user-context
1118000e7ae3SMatthew Knepley 
1119000e7ae3SMatthew Knepley   Level: developer
1120000e7ae3SMatthew Knepley 
1121000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1122000e7ae3SMatthew Knepley @*/
112363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1124000e7ae3SMatthew Knepley {
1125000e7ae3SMatthew Knepley   PetscFunctionBegin;
1126000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1127000e7ae3SMatthew Knepley }
1128000e7ae3SMatthew Knepley 
1129e74ef692SMatthew Knepley #undef __FUNCT__
1130e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSystemMatrixBC"
1131ac226902SBarry Smith /*@C
1132000e7ae3SMatthew Knepley   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1133000e7ae3SMatthew Knepley   to the system matrix and preconditioner of each system.
1134000e7ae3SMatthew Knepley 
1135000e7ae3SMatthew Knepley   Collective on TS
1136000e7ae3SMatthew Knepley 
1137000e7ae3SMatthew Knepley   Input Parameters:
1138000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1139000e7ae3SMatthew Knepley - func - The function
1140000e7ae3SMatthew Knepley 
1141000e7ae3SMatthew Knepley   Calling sequence of func:
1142000e7ae3SMatthew Knepley . func (TS ts, Mat A, Mat B, void *ctx);
1143000e7ae3SMatthew Knepley 
1144000e7ae3SMatthew Knepley + A   - The current system matrix
1145000e7ae3SMatthew Knepley . B   - The current preconditioner
1146000e7ae3SMatthew Knepley - ctx - The user-context
1147000e7ae3SMatthew Knepley 
1148000e7ae3SMatthew Knepley   Level: intermediate
1149000e7ae3SMatthew Knepley 
1150000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1151000e7ae3SMatthew Knepley @*/
115263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSystemMatrixBC(TS ts, PetscErrorCode (*func)(TS, Mat, Mat, void *))
1153000e7ae3SMatthew Knepley {
1154000e7ae3SMatthew Knepley   PetscFunctionBegin;
11554482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1156000e7ae3SMatthew Knepley   ts->ops->applymatrixbc = func;
1157000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1158000e7ae3SMatthew Knepley }
1159000e7ae3SMatthew Knepley 
1160e74ef692SMatthew Knepley #undef __FUNCT__
1161e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSystemMatrixBC"
1162000e7ae3SMatthew Knepley /*@
1163000e7ae3SMatthew Knepley   TSDefaultSystemMatrixBC - The default boundary condition function which
1164000e7ae3SMatthew Knepley   does nothing.
1165000e7ae3SMatthew Knepley 
1166000e7ae3SMatthew Knepley   Collective on TS
1167000e7ae3SMatthew Knepley 
1168000e7ae3SMatthew Knepley   Input Parameters:
1169000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1170000e7ae3SMatthew Knepley . A   - The system matrix
1171000e7ae3SMatthew Knepley . B   - The preconditioner
1172000e7ae3SMatthew Knepley - ctx - The user-context
1173000e7ae3SMatthew Knepley 
1174000e7ae3SMatthew Knepley   Level: developer
1175000e7ae3SMatthew Knepley 
1176000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1177000e7ae3SMatthew Knepley @*/
117863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1179000e7ae3SMatthew Knepley {
1180000e7ae3SMatthew Knepley   PetscFunctionBegin;
1181000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1182000e7ae3SMatthew Knepley }
1183000e7ae3SMatthew Knepley 
1184e74ef692SMatthew Knepley #undef __FUNCT__
1185e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSolutionBC"
1186ac226902SBarry Smith /*@C
1187000e7ae3SMatthew Knepley   TSSetSolutionBC - Sets the function which applies boundary conditions
1188000e7ae3SMatthew Knepley   to the solution of each system. This is necessary in nonlinear systems
1189000e7ae3SMatthew Knepley   which time dependent boundary conditions.
1190000e7ae3SMatthew Knepley 
1191000e7ae3SMatthew Knepley   Collective on TS
1192000e7ae3SMatthew Knepley 
1193000e7ae3SMatthew Knepley   Input Parameters:
1194000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1195000e7ae3SMatthew Knepley - func - The function
1196000e7ae3SMatthew Knepley 
1197000e7ae3SMatthew Knepley   Calling sequence of func:
1198000e7ae3SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
1199000e7ae3SMatthew Knepley 
1200000e7ae3SMatthew Knepley + sol - The current solution vector
1201000e7ae3SMatthew Knepley - ctx - The user-context
1202000e7ae3SMatthew Knepley 
1203000e7ae3SMatthew Knepley   Level: intermediate
1204000e7ae3SMatthew Knepley 
1205000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1206000e7ae3SMatthew Knepley @*/
120763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSolutionBC(TS ts, PetscErrorCode (*func)(TS, Vec, void *))
1208000e7ae3SMatthew Knepley {
1209000e7ae3SMatthew Knepley   PetscFunctionBegin;
12104482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1211000e7ae3SMatthew Knepley   ts->ops->applysolbc = func;
1212000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1213000e7ae3SMatthew Knepley }
1214000e7ae3SMatthew Knepley 
1215e74ef692SMatthew Knepley #undef __FUNCT__
1216e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSolutionBC"
1217000e7ae3SMatthew Knepley /*@
1218000e7ae3SMatthew Knepley   TSDefaultSolutionBC - The default boundary condition function which
1219000e7ae3SMatthew Knepley   does nothing.
1220000e7ae3SMatthew Knepley 
1221000e7ae3SMatthew Knepley   Collective on TS
1222000e7ae3SMatthew Knepley 
1223000e7ae3SMatthew Knepley   Input Parameters:
1224000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1225000e7ae3SMatthew Knepley . sol - The solution
1226000e7ae3SMatthew Knepley - ctx - The user-context
1227000e7ae3SMatthew Knepley 
1228000e7ae3SMatthew Knepley   Level: developer
1229000e7ae3SMatthew Knepley 
1230000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1231000e7ae3SMatthew Knepley @*/
123263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1233000e7ae3SMatthew Knepley {
1234000e7ae3SMatthew Knepley   PetscFunctionBegin;
1235000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1236000e7ae3SMatthew Knepley }
1237000e7ae3SMatthew Knepley 
1238e74ef692SMatthew Knepley #undef __FUNCT__
1239e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1240ac226902SBarry Smith /*@C
1241000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1242000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1243000e7ae3SMatthew Knepley 
1244000e7ae3SMatthew Knepley   Collective on TS
1245000e7ae3SMatthew Knepley 
1246000e7ae3SMatthew Knepley   Input Parameters:
1247000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1248000e7ae3SMatthew Knepley - func - The function
1249000e7ae3SMatthew Knepley 
1250000e7ae3SMatthew Knepley   Calling sequence of func:
1251000e7ae3SMatthew Knepley . func (TS ts);
1252000e7ae3SMatthew Knepley 
1253000e7ae3SMatthew Knepley   Level: intermediate
1254000e7ae3SMatthew Knepley 
1255000e7ae3SMatthew Knepley .keywords: TS, timestep
1256000e7ae3SMatthew Knepley @*/
125763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1258000e7ae3SMatthew Knepley {
1259000e7ae3SMatthew Knepley   PetscFunctionBegin;
12604482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1261000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1262000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1263000e7ae3SMatthew Knepley }
1264000e7ae3SMatthew Knepley 
1265e74ef692SMatthew Knepley #undef __FUNCT__
1266e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1267000e7ae3SMatthew Knepley /*@
1268000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1269000e7ae3SMatthew Knepley 
1270000e7ae3SMatthew Knepley   Collective on TS
1271000e7ae3SMatthew Knepley 
1272000e7ae3SMatthew Knepley   Input Parameters:
1273000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1274000e7ae3SMatthew Knepley 
1275000e7ae3SMatthew Knepley   Level: developer
1276000e7ae3SMatthew Knepley 
1277000e7ae3SMatthew Knepley .keywords: TS, timestep
1278000e7ae3SMatthew Knepley @*/
127963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts)
1280000e7ae3SMatthew Knepley {
1281000e7ae3SMatthew Knepley   PetscFunctionBegin;
1282000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1283000e7ae3SMatthew Knepley }
1284000e7ae3SMatthew Knepley 
1285e74ef692SMatthew Knepley #undef __FUNCT__
1286e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1287ac226902SBarry Smith /*@C
1288000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1289000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1290000e7ae3SMatthew Knepley   the time step.
1291000e7ae3SMatthew Knepley 
1292000e7ae3SMatthew Knepley   Collective on TS
1293000e7ae3SMatthew Knepley 
1294000e7ae3SMatthew Knepley   Input Parameters:
1295000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1296000e7ae3SMatthew Knepley - func - The function
1297000e7ae3SMatthew Knepley 
1298000e7ae3SMatthew Knepley   Calling sequence of func:
1299000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1300000e7ae3SMatthew Knepley 
1301000e7ae3SMatthew Knepley + t   - The current time
1302000e7ae3SMatthew Knepley - dt  - The current time step
1303000e7ae3SMatthew Knepley 
1304000e7ae3SMatthew Knepley   Level: intermediate
1305000e7ae3SMatthew Knepley 
1306000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1307000e7ae3SMatthew Knepley @*/
130863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUpdate(TS ts, PetscErrorCode (*func)(TS, PetscReal, PetscReal *))
1309000e7ae3SMatthew Knepley {
1310000e7ae3SMatthew Knepley   PetscFunctionBegin;
13114482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1312000e7ae3SMatthew Knepley   ts->ops->update = func;
1313000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1314000e7ae3SMatthew Knepley }
1315000e7ae3SMatthew Knepley 
1316e74ef692SMatthew Knepley #undef __FUNCT__
1317e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1318000e7ae3SMatthew Knepley /*@
1319000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1320000e7ae3SMatthew Knepley 
1321000e7ae3SMatthew Knepley   Collective on TS
1322000e7ae3SMatthew Knepley 
1323000e7ae3SMatthew Knepley   Input Parameters:
1324000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1325000e7ae3SMatthew Knepley - t   - The current time
1326000e7ae3SMatthew Knepley 
1327000e7ae3SMatthew Knepley   Output Parameters:
1328000e7ae3SMatthew Knepley . dt  - The current time step
1329000e7ae3SMatthew Knepley 
1330000e7ae3SMatthew Knepley   Level: developer
1331000e7ae3SMatthew Knepley 
1332000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1333000e7ae3SMatthew Knepley @*/
133463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1335000e7ae3SMatthew Knepley {
1336000e7ae3SMatthew Knepley   PetscFunctionBegin;
1337000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1338000e7ae3SMatthew Knepley }
1339000e7ae3SMatthew Knepley 
1340e74ef692SMatthew Knepley #undef __FUNCT__
1341e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1342ac226902SBarry Smith /*@C
1343000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1344000e7ae3SMatthew Knepley   called once at the end of time stepping.
1345000e7ae3SMatthew Knepley 
1346000e7ae3SMatthew Knepley   Collective on TS
1347000e7ae3SMatthew Knepley 
1348000e7ae3SMatthew Knepley   Input Parameters:
1349000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1350000e7ae3SMatthew Knepley - func - The function
1351000e7ae3SMatthew Knepley 
1352000e7ae3SMatthew Knepley   Calling sequence of func:
1353000e7ae3SMatthew Knepley . func (TS ts);
1354000e7ae3SMatthew Knepley 
1355000e7ae3SMatthew Knepley   Level: intermediate
1356000e7ae3SMatthew Knepley 
1357000e7ae3SMatthew Knepley .keywords: TS, timestep
1358000e7ae3SMatthew Knepley @*/
135963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1360000e7ae3SMatthew Knepley {
1361000e7ae3SMatthew Knepley   PetscFunctionBegin;
13624482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1363000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1364000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1365000e7ae3SMatthew Knepley }
1366000e7ae3SMatthew Knepley 
1367e74ef692SMatthew Knepley #undef __FUNCT__
1368e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1369000e7ae3SMatthew Knepley /*@
1370000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1371000e7ae3SMatthew Knepley 
1372000e7ae3SMatthew Knepley   Collective on TS
1373000e7ae3SMatthew Knepley 
1374000e7ae3SMatthew Knepley   Input Parameters:
1375000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1376000e7ae3SMatthew Knepley 
1377000e7ae3SMatthew Knepley   Level: developer
1378000e7ae3SMatthew Knepley 
1379000e7ae3SMatthew Knepley .keywords: TS, timestep
1380000e7ae3SMatthew Knepley @*/
138163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts)
1382000e7ae3SMatthew Knepley {
1383000e7ae3SMatthew Knepley   PetscFunctionBegin;
1384000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1385000e7ae3SMatthew Knepley }
1386000e7ae3SMatthew Knepley 
1387d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1388d763cef2SBarry Smith 
13894a2ae208SSatish Balay #undef __FUNCT__
13904a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
1391d763cef2SBarry Smith /*@C
1392d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1393d763cef2SBarry Smith    timestep to display the iteration's  progress.
1394d763cef2SBarry Smith 
1395d763cef2SBarry Smith    Collective on TS
1396d763cef2SBarry Smith 
1397d763cef2SBarry Smith    Input Parameters:
1398d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1399d763cef2SBarry Smith .  func - monitoring routine
1400329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1401b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1402b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1403b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1404d763cef2SBarry Smith 
1405d763cef2SBarry Smith    Calling sequence of func:
1406a7cc72afSBarry Smith $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1407d763cef2SBarry Smith 
1408d763cef2SBarry Smith +    ts - the TS context
1409d763cef2SBarry Smith .    steps - iteration number
14101f06c33eSBarry Smith .    time - current time
1411d763cef2SBarry Smith .    x - current iterate
1412d763cef2SBarry Smith -    mctx - [optional] monitoring context
1413d763cef2SBarry Smith 
1414d763cef2SBarry Smith    Notes:
1415d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1416d763cef2SBarry Smith    already has been loaded.
1417d763cef2SBarry Smith 
1418d763cef2SBarry Smith    Level: intermediate
1419d763cef2SBarry Smith 
1420d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1421d763cef2SBarry Smith 
1422d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
1423d763cef2SBarry Smith @*/
142463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetMonitor(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1425d763cef2SBarry Smith {
1426d763cef2SBarry Smith   PetscFunctionBegin;
14274482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1428d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
142929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1430d763cef2SBarry Smith   }
1431d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1432329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1433d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1434d763cef2SBarry Smith   PetscFunctionReturn(0);
1435d763cef2SBarry Smith }
1436d763cef2SBarry Smith 
14374a2ae208SSatish Balay #undef __FUNCT__
14384a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
1439d763cef2SBarry Smith /*@C
1440d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1441d763cef2SBarry Smith 
1442d763cef2SBarry Smith    Collective on TS
1443d763cef2SBarry Smith 
1444d763cef2SBarry Smith    Input Parameters:
1445d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1446d763cef2SBarry Smith 
1447d763cef2SBarry Smith    Notes:
1448d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1449d763cef2SBarry Smith 
1450d763cef2SBarry Smith    Level: intermediate
1451d763cef2SBarry Smith 
1452d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1453d763cef2SBarry Smith 
1454d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
1455d763cef2SBarry Smith @*/
145663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSClearMonitor(TS ts)
1457d763cef2SBarry Smith {
1458d952e501SBarry Smith   PetscErrorCode ierr;
1459d952e501SBarry Smith   PetscInt       i;
1460d952e501SBarry Smith 
1461d763cef2SBarry Smith   PetscFunctionBegin;
14624482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1463d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1464d952e501SBarry Smith     if (ts->mdestroy[i]) {
1465d952e501SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
1466d952e501SBarry Smith     }
1467d952e501SBarry Smith   }
1468d763cef2SBarry Smith   ts->numbermonitors = 0;
1469d763cef2SBarry Smith   PetscFunctionReturn(0);
1470d763cef2SBarry Smith }
1471d763cef2SBarry Smith 
14724a2ae208SSatish Balay #undef __FUNCT__
14734a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
1474d8e5e3e6SSatish Balay /*@
1475d8e5e3e6SSatish Balay    TSDefaultMonitor - Sets the Default monitor
1476d8e5e3e6SSatish Balay @*/
1477a7cc72afSBarry Smith PetscErrorCode TSDefaultMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1478d763cef2SBarry Smith {
1479dfbe8321SBarry Smith   PetscErrorCode ierr;
1480eabae89aSBarry Smith   PetscViewer    viewer = (PetscViewer)ctx;
1481d132466eSBarry Smith 
1482d763cef2SBarry Smith   PetscFunctionBegin;
1483eabae89aSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
1484a83599f4SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"timestep %D dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1485d763cef2SBarry Smith   PetscFunctionReturn(0);
1486d763cef2SBarry Smith }
1487d763cef2SBarry Smith 
14884a2ae208SSatish Balay #undef __FUNCT__
14894a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1490d763cef2SBarry Smith /*@
1491d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1492d763cef2SBarry Smith 
1493d763cef2SBarry Smith    Collective on TS
1494d763cef2SBarry Smith 
1495d763cef2SBarry Smith    Input Parameter:
1496d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1497d763cef2SBarry Smith 
1498d763cef2SBarry Smith    Output Parameters:
1499d763cef2SBarry Smith +  steps - number of iterations until termination
1500142b95e3SSatish Balay -  ptime - time until termination
1501d763cef2SBarry Smith 
1502d763cef2SBarry Smith    Level: beginner
1503d763cef2SBarry Smith 
1504d763cef2SBarry Smith .keywords: TS, timestep, solve
1505d763cef2SBarry Smith 
1506d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1507d763cef2SBarry Smith @*/
150863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1509d763cef2SBarry Smith {
1510dfbe8321SBarry Smith   PetscErrorCode ierr;
1511d763cef2SBarry Smith 
1512d763cef2SBarry Smith   PetscFunctionBegin;
15134482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1514d405a339SMatthew Knepley   if (!ts->setupcalled) {
1515d405a339SMatthew Knepley     ierr = TSSetUp(ts);CHKERRQ(ierr);
1516d405a339SMatthew Knepley   }
1517d405a339SMatthew Knepley 
1518d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1519d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1520000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1521d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1522d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1523d405a339SMatthew Knepley 
15244bb05414SBarry Smith   if (!PetscPreLoadingOn) {
15254bb05414SBarry Smith     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1526d405a339SMatthew Knepley   }
1527d763cef2SBarry Smith   PetscFunctionReturn(0);
1528d763cef2SBarry Smith }
1529d763cef2SBarry Smith 
15304a2ae208SSatish Balay #undef __FUNCT__
15314a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1532d763cef2SBarry Smith /*
1533d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1534d763cef2SBarry Smith */
1535a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1536d763cef2SBarry Smith {
15376849ba73SBarry Smith   PetscErrorCode ierr;
1538a7cc72afSBarry Smith   PetscInt i,n = ts->numbermonitors;
1539d763cef2SBarry Smith 
1540d763cef2SBarry Smith   PetscFunctionBegin;
1541d763cef2SBarry Smith   for (i=0; i<n; i++) {
1542142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1543d763cef2SBarry Smith   }
1544d763cef2SBarry Smith   PetscFunctionReturn(0);
1545d763cef2SBarry Smith }
1546d763cef2SBarry Smith 
1547d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1548d763cef2SBarry Smith 
15494a2ae208SSatish Balay #undef __FUNCT__
15504a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1551d763cef2SBarry Smith /*@C
1552d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1553d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1554d763cef2SBarry Smith 
1555d763cef2SBarry Smith    Collective on TS
1556d763cef2SBarry Smith 
1557d763cef2SBarry Smith    Input Parameters:
1558d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1559d763cef2SBarry Smith .  label - the title to put in the title bar
15607c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1561d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1562d763cef2SBarry Smith 
1563d763cef2SBarry Smith    Output Parameter:
1564d763cef2SBarry Smith .  draw - the drawing context
1565d763cef2SBarry Smith 
1566d763cef2SBarry Smith    Options Database Key:
1567d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1568d763cef2SBarry Smith 
1569d763cef2SBarry Smith    Notes:
1570b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1571d763cef2SBarry Smith 
1572d763cef2SBarry Smith    Level: intermediate
1573d763cef2SBarry Smith 
15747c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1575d763cef2SBarry Smith 
1576d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
15777c922b88SBarry Smith 
1578d763cef2SBarry Smith @*/
157963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1580d763cef2SBarry Smith {
1581b0a32e0cSBarry Smith   PetscDraw      win;
1582dfbe8321SBarry Smith   PetscErrorCode ierr;
1583d763cef2SBarry Smith 
1584d763cef2SBarry Smith   PetscFunctionBegin;
1585b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1586b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1587b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1588b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1589d763cef2SBarry Smith 
159052e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1591d763cef2SBarry Smith   PetscFunctionReturn(0);
1592d763cef2SBarry Smith }
1593d763cef2SBarry Smith 
15944a2ae208SSatish Balay #undef __FUNCT__
15954a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
1596a7cc72afSBarry Smith PetscErrorCode TSLGMonitor(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1597d763cef2SBarry Smith {
1598b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
159987828ca2SBarry Smith   PetscReal      x,y = ptime;
1600dfbe8321SBarry Smith   PetscErrorCode ierr;
1601d763cef2SBarry Smith 
1602d763cef2SBarry Smith   PetscFunctionBegin;
16037c922b88SBarry Smith   if (!monctx) {
16047c922b88SBarry Smith     MPI_Comm    comm;
1605b0a32e0cSBarry Smith     PetscViewer viewer;
16067c922b88SBarry Smith 
16077c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1608b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1609b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
16107c922b88SBarry Smith   }
16117c922b88SBarry Smith 
1612b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
161387828ca2SBarry Smith   x = (PetscReal)n;
1614b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1615d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1616b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1617d763cef2SBarry Smith   }
1618d763cef2SBarry Smith   PetscFunctionReturn(0);
1619d763cef2SBarry Smith }
1620d763cef2SBarry Smith 
16214a2ae208SSatish Balay #undef __FUNCT__
16224a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1623d763cef2SBarry Smith /*@C
1624d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1625d763cef2SBarry Smith    with TSLGMonitorCreate().
1626d763cef2SBarry Smith 
1627b0a32e0cSBarry Smith    Collective on PetscDrawLG
1628d763cef2SBarry Smith 
1629d763cef2SBarry Smith    Input Parameter:
1630d763cef2SBarry Smith .  draw - the drawing context
1631d763cef2SBarry Smith 
1632d763cef2SBarry Smith    Level: intermediate
1633d763cef2SBarry Smith 
1634d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1635d763cef2SBarry Smith 
1636d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1637d763cef2SBarry Smith @*/
163863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSLGMonitorDestroy(PetscDrawLG drawlg)
1639d763cef2SBarry Smith {
1640b0a32e0cSBarry Smith   PetscDraw      draw;
1641dfbe8321SBarry Smith   PetscErrorCode ierr;
1642d763cef2SBarry Smith 
1643d763cef2SBarry Smith   PetscFunctionBegin;
1644b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1645b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1646b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1647d763cef2SBarry Smith   PetscFunctionReturn(0);
1648d763cef2SBarry Smith }
1649d763cef2SBarry Smith 
16504a2ae208SSatish Balay #undef __FUNCT__
16514a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1652d763cef2SBarry Smith /*@
1653d763cef2SBarry Smith    TSGetTime - Gets the current time.
1654d763cef2SBarry Smith 
1655d763cef2SBarry Smith    Not Collective
1656d763cef2SBarry Smith 
1657d763cef2SBarry Smith    Input Parameter:
1658d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1659d763cef2SBarry Smith 
1660d763cef2SBarry Smith    Output Parameter:
1661d763cef2SBarry Smith .  t  - the current time
1662d763cef2SBarry Smith 
1663d763cef2SBarry Smith    Contributed by: Matthew Knepley
1664d763cef2SBarry Smith 
1665d763cef2SBarry Smith    Level: beginner
1666d763cef2SBarry Smith 
1667d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1668d763cef2SBarry Smith 
1669d763cef2SBarry Smith .keywords: TS, get, time
1670d763cef2SBarry Smith @*/
167163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t)
1672d763cef2SBarry Smith {
1673d763cef2SBarry Smith   PetscFunctionBegin;
16744482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
16754482741eSBarry Smith   PetscValidDoublePointer(t,2);
1676d763cef2SBarry Smith   *t = ts->ptime;
1677d763cef2SBarry Smith   PetscFunctionReturn(0);
1678d763cef2SBarry Smith }
1679d763cef2SBarry Smith 
16804a2ae208SSatish Balay #undef __FUNCT__
16814a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1682d763cef2SBarry Smith /*@C
1683d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1684d763cef2SBarry Smith    TS options in the database.
1685d763cef2SBarry Smith 
1686d763cef2SBarry Smith    Collective on TS
1687d763cef2SBarry Smith 
1688d763cef2SBarry Smith    Input Parameter:
1689d763cef2SBarry Smith +  ts     - The TS context
1690d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1691d763cef2SBarry Smith 
1692d763cef2SBarry Smith    Notes:
1693d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1694d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1695d763cef2SBarry Smith    hyphen.
1696d763cef2SBarry Smith 
1697d763cef2SBarry Smith    Contributed by: Matthew Knepley
1698d763cef2SBarry Smith 
1699d763cef2SBarry Smith    Level: advanced
1700d763cef2SBarry Smith 
1701d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1702d763cef2SBarry Smith 
1703d763cef2SBarry Smith .seealso: TSSetFromOptions()
1704d763cef2SBarry Smith 
1705d763cef2SBarry Smith @*/
170663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[])
1707d763cef2SBarry Smith {
1708dfbe8321SBarry Smith   PetscErrorCode ierr;
1709d763cef2SBarry Smith 
1710d763cef2SBarry Smith   PetscFunctionBegin;
17114482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1712d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1713d763cef2SBarry Smith   switch(ts->problem_type) {
1714d763cef2SBarry Smith     case TS_NONLINEAR:
1715d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1716d763cef2SBarry Smith       break;
1717d763cef2SBarry Smith     case TS_LINEAR:
171894b7f48cSBarry Smith       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1719d763cef2SBarry Smith       break;
1720d763cef2SBarry Smith   }
1721d763cef2SBarry Smith   PetscFunctionReturn(0);
1722d763cef2SBarry Smith }
1723d763cef2SBarry Smith 
1724d763cef2SBarry Smith 
17254a2ae208SSatish Balay #undef __FUNCT__
17264a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1727d763cef2SBarry Smith /*@C
1728d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1729d763cef2SBarry Smith    TS options in the database.
1730d763cef2SBarry Smith 
1731d763cef2SBarry Smith    Collective on TS
1732d763cef2SBarry Smith 
1733d763cef2SBarry Smith    Input Parameter:
1734d763cef2SBarry Smith +  ts     - The TS context
1735d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1736d763cef2SBarry Smith 
1737d763cef2SBarry Smith    Notes:
1738d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1739d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1740d763cef2SBarry Smith    hyphen.
1741d763cef2SBarry Smith 
1742d763cef2SBarry Smith    Contributed by: Matthew Knepley
1743d763cef2SBarry Smith 
1744d763cef2SBarry Smith    Level: advanced
1745d763cef2SBarry Smith 
1746d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1747d763cef2SBarry Smith 
1748d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1749d763cef2SBarry Smith 
1750d763cef2SBarry Smith @*/
175163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[])
1752d763cef2SBarry Smith {
1753dfbe8321SBarry Smith   PetscErrorCode ierr;
1754d763cef2SBarry Smith 
1755d763cef2SBarry Smith   PetscFunctionBegin;
17564482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1757d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1758d763cef2SBarry Smith   switch(ts->problem_type) {
1759d763cef2SBarry Smith     case TS_NONLINEAR:
1760d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1761d763cef2SBarry Smith       break;
1762d763cef2SBarry Smith     case TS_LINEAR:
176394b7f48cSBarry Smith       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1764d763cef2SBarry Smith       break;
1765d763cef2SBarry Smith   }
1766d763cef2SBarry Smith   PetscFunctionReturn(0);
1767d763cef2SBarry Smith }
1768d763cef2SBarry Smith 
17694a2ae208SSatish Balay #undef __FUNCT__
17704a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1771d763cef2SBarry Smith /*@C
1772d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1773d763cef2SBarry Smith    TS options in the database.
1774d763cef2SBarry Smith 
1775d763cef2SBarry Smith    Not Collective
1776d763cef2SBarry Smith 
1777d763cef2SBarry Smith    Input Parameter:
1778d763cef2SBarry Smith .  ts - The TS context
1779d763cef2SBarry Smith 
1780d763cef2SBarry Smith    Output Parameter:
1781d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1782d763cef2SBarry Smith 
1783d763cef2SBarry Smith    Contributed by: Matthew Knepley
1784d763cef2SBarry Smith 
1785d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1786d763cef2SBarry Smith    sufficient length to hold the prefix.
1787d763cef2SBarry Smith 
1788d763cef2SBarry Smith    Level: intermediate
1789d763cef2SBarry Smith 
1790d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1791d763cef2SBarry Smith 
1792d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1793d763cef2SBarry Smith @*/
1794e060cb09SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[])
1795d763cef2SBarry Smith {
1796dfbe8321SBarry Smith   PetscErrorCode ierr;
1797d763cef2SBarry Smith 
1798d763cef2SBarry Smith   PetscFunctionBegin;
17994482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
18004482741eSBarry Smith   PetscValidPointer(prefix,2);
1801d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1802d763cef2SBarry Smith   PetscFunctionReturn(0);
1803d763cef2SBarry Smith }
1804d763cef2SBarry Smith 
18054a2ae208SSatish Balay #undef __FUNCT__
18064a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1807d763cef2SBarry Smith /*@C
1808d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1809d763cef2SBarry Smith 
1810d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1811d763cef2SBarry Smith 
1812d763cef2SBarry Smith    Input Parameter:
1813d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1814d763cef2SBarry Smith 
1815d763cef2SBarry Smith    Output Parameters:
1816d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1817d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1818d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1819d763cef2SBarry Smith 
1820d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1821d763cef2SBarry Smith 
1822d763cef2SBarry Smith    Contributed by: Matthew Knepley
1823d763cef2SBarry Smith 
1824d763cef2SBarry Smith    Level: intermediate
1825d763cef2SBarry Smith 
1826d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1827d763cef2SBarry Smith 
1828d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1829d763cef2SBarry Smith 
1830d763cef2SBarry Smith @*/
183163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1832d763cef2SBarry Smith {
1833d763cef2SBarry Smith   PetscFunctionBegin;
18344482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1835d763cef2SBarry Smith   if (A)   *A = ts->A;
1836d763cef2SBarry Smith   if (M)   *M = ts->B;
1837d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1838d763cef2SBarry Smith   PetscFunctionReturn(0);
1839d763cef2SBarry Smith }
1840d763cef2SBarry Smith 
18414a2ae208SSatish Balay #undef __FUNCT__
18424a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1843d763cef2SBarry Smith /*@C
1844d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1845d763cef2SBarry Smith 
1846d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1847d763cef2SBarry Smith 
1848d763cef2SBarry Smith    Input Parameter:
1849d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1850d763cef2SBarry Smith 
1851d763cef2SBarry Smith    Output Parameters:
1852d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1853d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1854d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1855d763cef2SBarry Smith 
1856d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1857d763cef2SBarry Smith 
1858d763cef2SBarry Smith    Contributed by: Matthew Knepley
1859d763cef2SBarry Smith 
1860d763cef2SBarry Smith    Level: intermediate
1861d763cef2SBarry Smith 
1862d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1863d763cef2SBarry Smith 
1864d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1865d763cef2SBarry Smith @*/
186663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1867d763cef2SBarry Smith {
1868dfbe8321SBarry Smith   PetscErrorCode ierr;
1869d763cef2SBarry Smith 
1870d763cef2SBarry Smith   PetscFunctionBegin;
1871d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1872d763cef2SBarry Smith   PetscFunctionReturn(0);
1873d763cef2SBarry Smith }
1874d763cef2SBarry Smith 
18751713a123SBarry Smith #undef __FUNCT__
18761713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
18771713a123SBarry Smith /*@C
18781713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
18791713a123SBarry Smith    VecView() for the solution at each timestep
18801713a123SBarry Smith 
18811713a123SBarry Smith    Collective on TS
18821713a123SBarry Smith 
18831713a123SBarry Smith    Input Parameters:
18841713a123SBarry Smith +  ts - the TS context
18851713a123SBarry Smith .  step - current time-step
1886142b95e3SSatish Balay .  ptime - current time
18871713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
18881713a123SBarry Smith 
18891713a123SBarry Smith    Level: intermediate
18901713a123SBarry Smith 
18911713a123SBarry Smith .keywords: TS,  vector, monitor, view
18921713a123SBarry Smith 
18931713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
18941713a123SBarry Smith @*/
189563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSVecViewMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
18961713a123SBarry Smith {
1897dfbe8321SBarry Smith   PetscErrorCode ierr;
18981713a123SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
18991713a123SBarry Smith 
19001713a123SBarry Smith   PetscFunctionBegin;
19011713a123SBarry Smith   if (!viewer) {
19021713a123SBarry Smith     MPI_Comm comm;
19031713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
19041713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
19051713a123SBarry Smith   }
19061713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
19071713a123SBarry Smith   PetscFunctionReturn(0);
19081713a123SBarry Smith }
19091713a123SBarry Smith 
19101713a123SBarry Smith 
19111713a123SBarry Smith 
1912