xref: /petsc/src/ts/interface/ts.c (revision a79aaaed103246988dafa581a233ebf306fec680)
1e090d566SSatish Balay #include "src/ts/tsimpl.h"        /*I "petscts.h"  I*/
2d763cef2SBarry Smith 
3d5ba7fb7SMatthew Knepley /* Logging support */
46849ba73SBarry Smith PetscCookie TS_COOKIE = 0;
56849ba73SBarry Smith PetscEvent  TS_Step = 0, TS_PseudoComputeTimeStep = 0, TS_FunctionEval = 0, TS_JacobianEval = 0;
6d405a339SMatthew Knepley 
74a2ae208SSatish Balay #undef __FUNCT__
8bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
9bdad233fSMatthew Knepley /*
10bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
11bdad233fSMatthew Knepley 
12bdad233fSMatthew Knepley   Collective on TS
13bdad233fSMatthew Knepley 
14bdad233fSMatthew Knepley   Input Parameter:
15bdad233fSMatthew Knepley . ts - The ts
16bdad233fSMatthew Knepley 
17bdad233fSMatthew Knepley   Level: intermediate
18bdad233fSMatthew Knepley 
19bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
20bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
21bdad233fSMatthew Knepley */
226849ba73SBarry Smith static PetscErrorCode TSSetTypeFromOptions(TS ts)
23bdad233fSMatthew Knepley {
24bdad233fSMatthew Knepley   PetscTruth     opt;
252fc52814SBarry Smith   const char     *defaultType;
26bdad233fSMatthew Knepley   char           typeName[256];
27dfbe8321SBarry Smith   PetscErrorCode ierr;
28bdad233fSMatthew Knepley 
29bdad233fSMatthew Knepley   PetscFunctionBegin;
30abc0a331SBarry Smith   if (ts->type_name) {
31bdad233fSMatthew Knepley     defaultType = ts->type_name;
32bdad233fSMatthew Knepley   } else {
33bdad233fSMatthew Knepley     defaultType = TS_EULER;
34bdad233fSMatthew Knepley   }
35bdad233fSMatthew Knepley 
36bdad233fSMatthew Knepley   if (!TSRegisterAllCalled) {
37bdad233fSMatthew Knepley     ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);
38bdad233fSMatthew Knepley   }
39bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
40a7cc72afSBarry Smith   if (opt) {
41bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
42bdad233fSMatthew Knepley   } else {
43bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
44bdad233fSMatthew Knepley   }
45bdad233fSMatthew Knepley   PetscFunctionReturn(0);
46bdad233fSMatthew Knepley }
47bdad233fSMatthew Knepley 
48bdad233fSMatthew Knepley #undef __FUNCT__
49bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
50bdad233fSMatthew Knepley /*@
51bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
52bdad233fSMatthew Knepley 
53bdad233fSMatthew Knepley    Collective on TS
54bdad233fSMatthew Knepley 
55bdad233fSMatthew Knepley    Input Parameter:
56bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
57bdad233fSMatthew Knepley 
58bdad233fSMatthew Knepley    Options Database Keys:
59bdad233fSMatthew Knepley +  -ts_type <type> - TS_EULER, TS_BEULER, TS_PVODE, TS_PSEUDO, TS_CRANK_NICHOLSON
60bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
61bdad233fSMatthew Knepley .  -ts_max_time time - maximum time to compute to
62bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
63bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
64bdad233fSMatthew Knepley -  -ts_xmonitor - plot information at each timestep
65bdad233fSMatthew Knepley 
66bdad233fSMatthew Knepley    Level: beginner
67bdad233fSMatthew Knepley 
68bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
69bdad233fSMatthew Knepley 
70bdad233fSMatthew Knepley .seealso: TSGetType
71bdad233fSMatthew Knepley @*/
72dfbe8321SBarry Smith PetscErrorCode TSSetFromOptions(TS ts)
73bdad233fSMatthew Knepley {
74bdad233fSMatthew Knepley   PetscReal      dt;
75bdad233fSMatthew Knepley   PetscTruth     opt;
76dfbe8321SBarry Smith   PetscErrorCode ierr;
77bdad233fSMatthew Knepley 
78bdad233fSMatthew Knepley   PetscFunctionBegin;
794482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
80bdad233fSMatthew Knepley   ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");CHKERRQ(ierr);
81bdad233fSMatthew Knepley 
82bdad233fSMatthew Knepley   /* Handle generic TS options */
83bdad233fSMatthew Knepley   ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
84bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
85bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
86bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
87a7cc72afSBarry Smith   if (opt) {
88bdad233fSMatthew Knepley     ts->initial_time_step = ts->time_step = dt;
89bdad233fSMatthew Knepley   }
90bdad233fSMatthew Knepley 
91bdad233fSMatthew Knepley   /* Monitor options */
92bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_monitor","Monitor timestep size","TSDefaultMonitor",&opt);CHKERRQ(ierr);
93a7cc72afSBarry Smith     if (opt) {
94bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
95bdad233fSMatthew Knepley     }
96bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_xmonitor","Monitor timestep size graphically","TSLGMonitor",&opt);CHKERRQ(ierr);
97a7cc72afSBarry Smith     if (opt) {
98bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
99bdad233fSMatthew Knepley     }
100bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_vecmonitor","Monitor solution graphically","TSVecViewMonitor",&opt);CHKERRQ(ierr);
101a7cc72afSBarry Smith     if (opt) {
102bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSVecViewMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
103bdad233fSMatthew Knepley     }
104bdad233fSMatthew Knepley 
105bdad233fSMatthew Knepley   /* Handle TS type options */
106bdad233fSMatthew Knepley   ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
107bdad233fSMatthew Knepley 
108bdad233fSMatthew Knepley   /* Handle specific TS options */
109abc0a331SBarry Smith   if (ts->ops->setfromoptions) {
110bdad233fSMatthew Knepley     ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
111bdad233fSMatthew Knepley   }
112bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
113bdad233fSMatthew Knepley 
114bdad233fSMatthew Knepley   /* Handle subobject options */
115bdad233fSMatthew Knepley   switch(ts->problem_type) {
116156fc9a6SMatthew Knepley     /* Should check for implicit/explicit */
117bdad233fSMatthew Knepley   case TS_LINEAR:
118abc0a331SBarry Smith     if (ts->ksp) {
1197c236d22SBarry Smith       ierr = KSPSetOperators(ts->ksp,ts->A,ts->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
12094b7f48cSBarry Smith       ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr);
121156fc9a6SMatthew Knepley     }
122bdad233fSMatthew Knepley     break;
123bdad233fSMatthew Knepley   case TS_NONLINEAR:
124abc0a331SBarry Smith     if (ts->snes) {
1257c236d22SBarry Smith       /* this is a bit of a hack, but it gets the matrix information into SNES earlier
1267c236d22SBarry Smith          so that SNES and KSP have more information to pick reasonable defaults
1277c236d22SBarry Smith          before they allow users to set options */
1287c236d22SBarry Smith       ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,0,ts);CHKERRQ(ierr);
129bdad233fSMatthew Knepley       ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
130156fc9a6SMatthew Knepley     }
131bdad233fSMatthew Knepley     break;
132bdad233fSMatthew Knepley   default:
13377431f27SBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type);
134bdad233fSMatthew Knepley   }
135bdad233fSMatthew Knepley 
136bdad233fSMatthew Knepley   PetscFunctionReturn(0);
137bdad233fSMatthew Knepley }
138bdad233fSMatthew Knepley 
139bdad233fSMatthew Knepley #undef  __FUNCT__
140bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
141bdad233fSMatthew Knepley /*@
142bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
143bdad233fSMatthew Knepley 
144bdad233fSMatthew Knepley   Collective on TS
145bdad233fSMatthew Knepley 
146bdad233fSMatthew Knepley   Input Parameter:
147bdad233fSMatthew Knepley . ts - The ts
148bdad233fSMatthew Knepley 
149bdad233fSMatthew Knepley   Level: intermediate
150bdad233fSMatthew Knepley 
151bdad233fSMatthew Knepley .keywords: TS, view, options, database
152bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
153bdad233fSMatthew Knepley @*/
154dfbe8321SBarry Smith PetscErrorCode TSViewFromOptions(TS ts,const char title[])
155bdad233fSMatthew Knepley {
156bdad233fSMatthew Knepley   PetscViewer    viewer;
157bdad233fSMatthew Knepley   PetscDraw      draw;
158bdad233fSMatthew Knepley   PetscTruth     opt;
159bdad233fSMatthew Knepley   char           typeName[1024];
160e10c49a3SBarry Smith   char           fileName[PETSC_MAX_PATH_LEN];
161de4209c5SBarry Smith   size_t         len;
162dfbe8321SBarry Smith   PetscErrorCode ierr;
163bdad233fSMatthew Knepley 
164bdad233fSMatthew Knepley   PetscFunctionBegin;
165bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);CHKERRQ(ierr);
166a7cc72afSBarry Smith   if (opt) {
167bdad233fSMatthew Knepley     ierr = PetscOptionsGetString(ts->prefix, "-ts_view", typeName, 1024, &opt);CHKERRQ(ierr);
168bdad233fSMatthew Knepley     ierr = PetscStrlen(typeName, &len);CHKERRQ(ierr);
169bdad233fSMatthew Knepley     if (len > 0) {
170bdad233fSMatthew Knepley       ierr = PetscViewerCreate(ts->comm, &viewer);CHKERRQ(ierr);
171bdad233fSMatthew Knepley       ierr = PetscViewerSetType(viewer, typeName);CHKERRQ(ierr);
172bdad233fSMatthew Knepley       ierr = PetscOptionsGetString(ts->prefix, "-ts_view_file", fileName, 1024, &opt);CHKERRQ(ierr);
173a7cc72afSBarry Smith       if (opt) {
174bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, fileName);CHKERRQ(ierr);
175bdad233fSMatthew Knepley       } else {
176bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, ts->name);CHKERRQ(ierr);
177bdad233fSMatthew Knepley       }
178bdad233fSMatthew Knepley       ierr = TSView(ts, viewer);CHKERRQ(ierr);
179bdad233fSMatthew Knepley       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
180bdad233fSMatthew Knepley       ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
181bdad233fSMatthew Knepley     } else {
182bdad233fSMatthew Knepley       ierr = TSView(ts, PETSC_NULL);CHKERRQ(ierr);
183bdad233fSMatthew Knepley     }
184bdad233fSMatthew Knepley   }
185bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);CHKERRQ(ierr);
186a7cc72afSBarry Smith   if (opt) {
187bdad233fSMatthew Knepley     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
188bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
189a7cc72afSBarry Smith     if (title) {
1901836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
191bdad233fSMatthew Knepley     } else {
192bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject) ts);CHKERRQ(ierr);
1931836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, ts->name);CHKERRQ(ierr);
194bdad233fSMatthew Knepley     }
195bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
196bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
197bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
198bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
199bdad233fSMatthew Knepley   }
200bdad233fSMatthew Knepley   PetscFunctionReturn(0);
201bdad233fSMatthew Knepley }
202bdad233fSMatthew Knepley 
203bdad233fSMatthew Knepley #undef __FUNCT__
2044a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
205a7a1495cSBarry Smith /*@
2068c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
207a7a1495cSBarry Smith       set with TSSetRHSJacobian().
208a7a1495cSBarry Smith 
209a7a1495cSBarry Smith    Collective on TS and Vec
210a7a1495cSBarry Smith 
211a7a1495cSBarry Smith    Input Parameters:
212a7a1495cSBarry Smith +  ts - the SNES context
213a7a1495cSBarry Smith .  t - current timestep
214a7a1495cSBarry Smith -  x - input vector
215a7a1495cSBarry Smith 
216a7a1495cSBarry Smith    Output Parameters:
217a7a1495cSBarry Smith +  A - Jacobian matrix
218a7a1495cSBarry Smith .  B - optional preconditioning matrix
219a7a1495cSBarry Smith -  flag - flag indicating matrix structure
220a7a1495cSBarry Smith 
221a7a1495cSBarry Smith    Notes:
222a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
223a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
224a7a1495cSBarry Smith 
22594b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
226a7a1495cSBarry Smith    flag parameter.
227a7a1495cSBarry Smith 
228a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
229a7a1495cSBarry Smith 
230a7a1495cSBarry Smith    Level: developer
231a7a1495cSBarry Smith 
232a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
233a7a1495cSBarry Smith 
23494b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
235a7a1495cSBarry Smith @*/
236dfbe8321SBarry Smith PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
237a7a1495cSBarry Smith {
238dfbe8321SBarry Smith   PetscErrorCode ierr;
239a7a1495cSBarry Smith 
240a7a1495cSBarry Smith   PetscFunctionBegin;
2414482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2424482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
243c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
244a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
24529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
246a7a1495cSBarry Smith   }
247000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
248d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
249a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
250a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
251000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
252a7a1495cSBarry Smith     PetscStackPop;
253d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
254a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2554482741eSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE,4);
2564482741eSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE,5);
257ef66eb69SBarry Smith   } else {
258ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260ef66eb69SBarry Smith     if (*A != *B) {
261ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
262ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
263ef66eb69SBarry Smith     }
264ef66eb69SBarry Smith   }
265a7a1495cSBarry Smith   PetscFunctionReturn(0);
266a7a1495cSBarry Smith }
267a7a1495cSBarry Smith 
2684a2ae208SSatish Balay #undef __FUNCT__
2694a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
270d763cef2SBarry Smith /*
271d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
272d763cef2SBarry Smith 
273d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
274d763cef2SBarry Smith    this routine applies the matrix.
275d763cef2SBarry Smith */
276dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
277d763cef2SBarry Smith {
278dfbe8321SBarry Smith   PetscErrorCode ierr;
279d763cef2SBarry Smith 
280d763cef2SBarry Smith   PetscFunctionBegin;
2814482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2824482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
2834482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
284d763cef2SBarry Smith 
285d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
286000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
287d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
288000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
289d763cef2SBarry Smith     PetscStackPop;
2907c922b88SBarry Smith   } else {
291000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
292d763cef2SBarry Smith       MatStructure flg;
293d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
294000e7ae3SMatthew Knepley       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
295d763cef2SBarry Smith       PetscStackPop;
296d763cef2SBarry Smith     }
297d763cef2SBarry Smith     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
2987c922b88SBarry Smith   }
299d763cef2SBarry Smith 
300d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
301d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
302d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
303d763cef2SBarry Smith 
304d763cef2SBarry Smith   PetscFunctionReturn(0);
305d763cef2SBarry Smith }
306d763cef2SBarry Smith 
3074a2ae208SSatish Balay #undef __FUNCT__
3084a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
309d763cef2SBarry Smith /*@C
310d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
311d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
312d763cef2SBarry Smith 
313d763cef2SBarry Smith     Collective on TS
314d763cef2SBarry Smith 
315d763cef2SBarry Smith     Input Parameters:
316d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
317d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
318d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
319d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
320d763cef2SBarry Smith 
321d763cef2SBarry Smith     Calling sequence of func:
32287828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
323d763cef2SBarry Smith 
324d763cef2SBarry Smith +   t - current timestep
325d763cef2SBarry Smith .   u - input vector
326d763cef2SBarry Smith .   F - function vector
327d763cef2SBarry Smith -   ctx - [optional] user-defined function context
328d763cef2SBarry Smith 
329d763cef2SBarry Smith     Important:
330d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
331d763cef2SBarry Smith 
332d763cef2SBarry Smith     Level: beginner
333d763cef2SBarry Smith 
334d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
335d763cef2SBarry Smith 
336d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
337d763cef2SBarry Smith @*/
3386849ba73SBarry Smith PetscErrorCode TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
339d763cef2SBarry Smith {
340d763cef2SBarry Smith   PetscFunctionBegin;
341d763cef2SBarry Smith 
3424482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
343d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
34429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
345d763cef2SBarry Smith   }
346000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
347d763cef2SBarry Smith   ts->funP             = ctx;
348d763cef2SBarry Smith   PetscFunctionReturn(0);
349d763cef2SBarry Smith }
350d763cef2SBarry Smith 
3514a2ae208SSatish Balay #undef __FUNCT__
3524a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
353d763cef2SBarry Smith /*@C
354d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
355d763cef2SBarry Smith    Also sets the location to store A.
356d763cef2SBarry Smith 
357d763cef2SBarry Smith    Collective on TS
358d763cef2SBarry Smith 
359d763cef2SBarry Smith    Input Parameters:
360d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
361d763cef2SBarry Smith .  A   - matrix
362d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
363d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
364d763cef2SBarry Smith          if A is not a function of t.
365d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
366d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
367d763cef2SBarry Smith 
368d763cef2SBarry Smith    Calling sequence of func:
369a7cc72afSBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
370d763cef2SBarry Smith 
371d763cef2SBarry Smith +  t - current timestep
372d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
373d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
374d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
37594b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
376d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
377d763cef2SBarry Smith 
378d763cef2SBarry Smith    Notes:
37994b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
380d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
381d763cef2SBarry Smith 
382d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
383d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
384d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
385d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
386d763cef2SBarry Smith    throughout the global iterations.
387d763cef2SBarry Smith 
388d763cef2SBarry Smith    Important:
389d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
390d763cef2SBarry Smith 
391d763cef2SBarry Smith    Level: beginner
392d763cef2SBarry Smith 
393d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
394d763cef2SBarry Smith 
395d763cef2SBarry Smith .seealso: TSSetRHSFunction()
396d763cef2SBarry Smith @*/
3976849ba73SBarry Smith PetscErrorCode TSSetRHSMatrix(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
398d763cef2SBarry Smith {
399d763cef2SBarry Smith   PetscFunctionBegin;
4004482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4014482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
4024482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
403c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
404c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
405d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
40629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
407d763cef2SBarry Smith   }
408d763cef2SBarry Smith 
409000e7ae3SMatthew Knepley   ts->ops->rhsmatrix = f;
410d763cef2SBarry Smith   ts->jacP           = ctx;
411d763cef2SBarry Smith   ts->A              = A;
412d763cef2SBarry Smith   ts->B              = B;
413d763cef2SBarry Smith 
414d763cef2SBarry Smith   PetscFunctionReturn(0);
415d763cef2SBarry Smith }
416d763cef2SBarry Smith 
4174a2ae208SSatish Balay #undef __FUNCT__
4184a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
419d763cef2SBarry Smith /*@C
420d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
421d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
4223c94ec11SBarry Smith    Use TSSetRHSMatrix() for linear problems.
423d763cef2SBarry Smith 
424d763cef2SBarry Smith    Collective on TS
425d763cef2SBarry Smith 
426d763cef2SBarry Smith    Input Parameters:
427d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
428d763cef2SBarry Smith .  A   - Jacobian matrix
429d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
430d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
431d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
432d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
433d763cef2SBarry Smith 
434d763cef2SBarry Smith    Calling sequence of func:
43587828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
436d763cef2SBarry Smith 
437d763cef2SBarry Smith +  t - current timestep
438d763cef2SBarry Smith .  u - input vector
439d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
440d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
441d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
44294b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
443d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
444d763cef2SBarry Smith 
445d763cef2SBarry Smith    Notes:
44694b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
447d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
448d763cef2SBarry Smith 
449d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
450d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
451d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
452d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
453d763cef2SBarry Smith    throughout the global iterations.
454d763cef2SBarry Smith 
455d763cef2SBarry Smith    Level: beginner
456d763cef2SBarry Smith 
457d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
458d763cef2SBarry Smith 
459d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
4603c94ec11SBarry Smith           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetRHSMatrix()
461d763cef2SBarry Smith 
462d763cef2SBarry Smith @*/
4636849ba73SBarry Smith PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
464d763cef2SBarry Smith {
465d763cef2SBarry Smith   PetscFunctionBegin;
4664482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4674482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
4684482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
469c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
470c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
471d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
47229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
473d763cef2SBarry Smith   }
474d763cef2SBarry Smith 
475000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
476d763cef2SBarry Smith   ts->jacP             = ctx;
477d763cef2SBarry Smith   ts->A                = A;
478d763cef2SBarry Smith   ts->B                = B;
479d763cef2SBarry Smith   PetscFunctionReturn(0);
480d763cef2SBarry Smith }
481d763cef2SBarry Smith 
4824a2ae208SSatish Balay #undef __FUNCT__
4834a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
484d763cef2SBarry Smith /*
485d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
486d763cef2SBarry Smith 
487d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
488d763cef2SBarry Smith    this routine applies the matrix.
489d763cef2SBarry Smith */
490dfbe8321SBarry Smith PetscErrorCode TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
491d763cef2SBarry Smith {
492dfbe8321SBarry Smith   PetscErrorCode ierr;
493d763cef2SBarry Smith 
494d763cef2SBarry Smith   PetscFunctionBegin;
4954482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4964482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
497c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,x,3);
498d763cef2SBarry Smith 
499000e7ae3SMatthew Knepley   if (ts->ops->rhsbc) {
500d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
501000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
502d763cef2SBarry Smith     PetscStackPop;
503d763cef2SBarry Smith     PetscFunctionReturn(0);
504d763cef2SBarry Smith   }
505d763cef2SBarry Smith 
506d763cef2SBarry Smith   PetscFunctionReturn(0);
507d763cef2SBarry Smith }
508d763cef2SBarry Smith 
5094a2ae208SSatish Balay #undef __FUNCT__
5104a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
511d763cef2SBarry Smith /*@C
512d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
513d763cef2SBarry Smith     boundary conditions for the function F.
514d763cef2SBarry Smith 
515d763cef2SBarry Smith     Collective on TS
516d763cef2SBarry Smith 
517d763cef2SBarry Smith     Input Parameters:
518d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
519d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
520d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
521d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
522d763cef2SBarry Smith 
523d763cef2SBarry Smith     Calling sequence of func:
52487828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
525d763cef2SBarry Smith 
526d763cef2SBarry Smith +   t - current timestep
527d763cef2SBarry Smith .   F - function vector
528d763cef2SBarry Smith -   ctx - [optional] user-defined function context
529d763cef2SBarry Smith 
530d763cef2SBarry Smith     Level: intermediate
531d763cef2SBarry Smith 
532d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
533d763cef2SBarry Smith @*/
5346849ba73SBarry Smith PetscErrorCode TSSetRHSBoundaryConditions(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
535d763cef2SBarry Smith {
536d763cef2SBarry Smith   PetscFunctionBegin;
537d763cef2SBarry Smith 
5384482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
539d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
54029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
541d763cef2SBarry Smith   }
542000e7ae3SMatthew Knepley   ts->ops->rhsbc = f;
543d763cef2SBarry Smith   ts->bcP        = ctx;
544d763cef2SBarry Smith   PetscFunctionReturn(0);
545d763cef2SBarry Smith }
546d763cef2SBarry Smith 
5474a2ae208SSatish Balay #undef __FUNCT__
5484a2ae208SSatish Balay #define __FUNCT__ "TSView"
5497e2c5f70SBarry Smith /*@C
550d763cef2SBarry Smith     TSView - Prints the TS data structure.
551d763cef2SBarry Smith 
5524c49b128SBarry Smith     Collective on TS
553d763cef2SBarry Smith 
554d763cef2SBarry Smith     Input Parameters:
555d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
556d763cef2SBarry Smith -   viewer - visualization context
557d763cef2SBarry Smith 
558d763cef2SBarry Smith     Options Database Key:
559d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
560d763cef2SBarry Smith 
561d763cef2SBarry Smith     Notes:
562d763cef2SBarry Smith     The available visualization contexts include
563b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
564b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
565d763cef2SBarry Smith          output where only the first processor opens
566d763cef2SBarry Smith          the file.  All other processors send their
567d763cef2SBarry Smith          data to the first processor to print.
568d763cef2SBarry Smith 
569d763cef2SBarry Smith     The user can open an alternative visualization context with
570b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
571d763cef2SBarry Smith 
572d763cef2SBarry Smith     Level: beginner
573d763cef2SBarry Smith 
574d763cef2SBarry Smith .keywords: TS, timestep, view
575d763cef2SBarry Smith 
576b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
577d763cef2SBarry Smith @*/
578dfbe8321SBarry Smith PetscErrorCode TSView(TS ts,PetscViewer viewer)
579d763cef2SBarry Smith {
580dfbe8321SBarry Smith   PetscErrorCode ierr;
581454a90a3SBarry Smith   char           *type;
58232077d6dSBarry Smith   PetscTruth     iascii,isstring;
583d763cef2SBarry Smith 
584d763cef2SBarry Smith   PetscFunctionBegin;
5854482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
586b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
5874482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
588c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
589fd16b177SBarry Smith 
59032077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
591b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
59232077d6dSBarry Smith   if (iascii) {
593b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
594454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
595454a90a3SBarry Smith     if (type) {
596b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
597184914b5SBarry Smith     } else {
598b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
599184914b5SBarry Smith     }
600000e7ae3SMatthew Knepley     if (ts->ops->view) {
601b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
602000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
603b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
604d763cef2SBarry Smith     }
60577431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
606b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
607d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
60877431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
609d763cef2SBarry Smith     }
61077431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
6110f5bd95cSBarry Smith   } else if (isstring) {
612454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
613b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
614d763cef2SBarry Smith   }
615b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
61694b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
617d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
618b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
619d763cef2SBarry Smith   PetscFunctionReturn(0);
620d763cef2SBarry Smith }
621d763cef2SBarry Smith 
622d763cef2SBarry Smith 
6234a2ae208SSatish Balay #undef __FUNCT__
6244a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
625d763cef2SBarry Smith /*@C
626d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
627d763cef2SBarry Smith    the timesteppers.
628d763cef2SBarry Smith 
629d763cef2SBarry Smith    Collective on TS
630d763cef2SBarry Smith 
631d763cef2SBarry Smith    Input Parameters:
632d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
633d763cef2SBarry Smith -  usrP - optional user context
634d763cef2SBarry Smith 
635d763cef2SBarry Smith    Level: intermediate
636d763cef2SBarry Smith 
637d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
638d763cef2SBarry Smith 
639d763cef2SBarry Smith .seealso: TSGetApplicationContext()
640d763cef2SBarry Smith @*/
641dfbe8321SBarry Smith PetscErrorCode TSSetApplicationContext(TS ts,void *usrP)
642d763cef2SBarry Smith {
643d763cef2SBarry Smith   PetscFunctionBegin;
6444482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
645d763cef2SBarry Smith   ts->user = usrP;
646d763cef2SBarry Smith   PetscFunctionReturn(0);
647d763cef2SBarry Smith }
648d763cef2SBarry Smith 
6494a2ae208SSatish Balay #undef __FUNCT__
6504a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
651d763cef2SBarry Smith /*@C
652d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
653d763cef2SBarry Smith     timestepper.
654d763cef2SBarry Smith 
655d763cef2SBarry Smith     Not Collective
656d763cef2SBarry Smith 
657d763cef2SBarry Smith     Input Parameter:
658d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
659d763cef2SBarry Smith 
660d763cef2SBarry Smith     Output Parameter:
661d763cef2SBarry Smith .   usrP - user context
662d763cef2SBarry Smith 
663d763cef2SBarry Smith     Level: intermediate
664d763cef2SBarry Smith 
665d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
666d763cef2SBarry Smith 
667d763cef2SBarry Smith .seealso: TSSetApplicationContext()
668d763cef2SBarry Smith @*/
669dfbe8321SBarry Smith PetscErrorCode TSGetApplicationContext(TS ts,void **usrP)
670d763cef2SBarry Smith {
671d763cef2SBarry Smith   PetscFunctionBegin;
6724482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
673d763cef2SBarry Smith   *usrP = ts->user;
674d763cef2SBarry Smith   PetscFunctionReturn(0);
675d763cef2SBarry Smith }
676d763cef2SBarry Smith 
6774a2ae208SSatish Balay #undef __FUNCT__
6784a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
679d763cef2SBarry Smith /*@
680d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
681d763cef2SBarry Smith 
682d763cef2SBarry Smith    Not Collective
683d763cef2SBarry Smith 
684d763cef2SBarry Smith    Input Parameter:
685d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
686d763cef2SBarry Smith 
687d763cef2SBarry Smith    Output Parameter:
688d763cef2SBarry Smith .  iter - number steps so far
689d763cef2SBarry Smith 
690d763cef2SBarry Smith    Level: intermediate
691d763cef2SBarry Smith 
692d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
693d763cef2SBarry Smith @*/
694a7cc72afSBarry Smith PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter)
695d763cef2SBarry Smith {
696d763cef2SBarry Smith   PetscFunctionBegin;
6974482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
6984482741eSBarry Smith   PetscValidIntPointer(iter,2);
699d763cef2SBarry Smith   *iter = ts->steps;
700d763cef2SBarry Smith   PetscFunctionReturn(0);
701d763cef2SBarry Smith }
702d763cef2SBarry Smith 
7034a2ae208SSatish Balay #undef __FUNCT__
7044a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
705d763cef2SBarry Smith /*@
706d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
707d763cef2SBarry Smith    as well as the initial time.
708d763cef2SBarry Smith 
709d763cef2SBarry Smith    Collective on TS
710d763cef2SBarry Smith 
711d763cef2SBarry Smith    Input Parameters:
712d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
713d763cef2SBarry Smith .  initial_time - the initial time
714d763cef2SBarry Smith -  time_step - the size of the timestep
715d763cef2SBarry Smith 
716d763cef2SBarry Smith    Level: intermediate
717d763cef2SBarry Smith 
718d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
719d763cef2SBarry Smith 
720d763cef2SBarry Smith .keywords: TS, set, initial, timestep
721d763cef2SBarry Smith @*/
722dfbe8321SBarry Smith PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
723d763cef2SBarry Smith {
724d763cef2SBarry Smith   PetscFunctionBegin;
7254482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
726d763cef2SBarry Smith   ts->time_step         = time_step;
727d763cef2SBarry Smith   ts->initial_time_step = time_step;
728d763cef2SBarry Smith   ts->ptime             = initial_time;
729d763cef2SBarry Smith   PetscFunctionReturn(0);
730d763cef2SBarry Smith }
731d763cef2SBarry Smith 
7324a2ae208SSatish Balay #undef __FUNCT__
7334a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
734d763cef2SBarry Smith /*@
735d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
736d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
737d763cef2SBarry Smith 
738d763cef2SBarry Smith    Collective on TS
739d763cef2SBarry Smith 
740d763cef2SBarry Smith    Input Parameters:
741d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
742d763cef2SBarry Smith -  time_step - the size of the timestep
743d763cef2SBarry Smith 
744d763cef2SBarry Smith    Level: intermediate
745d763cef2SBarry Smith 
746d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
747d763cef2SBarry Smith 
748d763cef2SBarry Smith .keywords: TS, set, timestep
749d763cef2SBarry Smith @*/
750dfbe8321SBarry Smith PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step)
751d763cef2SBarry Smith {
752d763cef2SBarry Smith   PetscFunctionBegin;
7534482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
754d763cef2SBarry Smith   ts->time_step = time_step;
755d763cef2SBarry Smith   PetscFunctionReturn(0);
756d763cef2SBarry Smith }
757d763cef2SBarry Smith 
7584a2ae208SSatish Balay #undef __FUNCT__
7594a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
760d763cef2SBarry Smith /*@
761d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
762d763cef2SBarry Smith 
763d763cef2SBarry Smith    Not Collective
764d763cef2SBarry Smith 
765d763cef2SBarry Smith    Input Parameter:
766d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
767d763cef2SBarry Smith 
768d763cef2SBarry Smith    Output Parameter:
769d763cef2SBarry Smith .  dt - the current timestep size
770d763cef2SBarry Smith 
771d763cef2SBarry Smith    Level: intermediate
772d763cef2SBarry Smith 
773d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
774d763cef2SBarry Smith 
775d763cef2SBarry Smith .keywords: TS, get, timestep
776d763cef2SBarry Smith @*/
777dfbe8321SBarry Smith PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt)
778d763cef2SBarry Smith {
779d763cef2SBarry Smith   PetscFunctionBegin;
7804482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
7814482741eSBarry Smith   PetscValidDoublePointer(dt,2);
782d763cef2SBarry Smith   *dt = ts->time_step;
783d763cef2SBarry Smith   PetscFunctionReturn(0);
784d763cef2SBarry Smith }
785d763cef2SBarry Smith 
7864a2ae208SSatish Balay #undef __FUNCT__
7874a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
788d763cef2SBarry Smith /*@C
789d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
790d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
791d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
792d763cef2SBarry Smith    the solution at the next timestep has been calculated.
793d763cef2SBarry Smith 
794d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
795d763cef2SBarry Smith 
796d763cef2SBarry Smith    Input Parameter:
797d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
798d763cef2SBarry Smith 
799d763cef2SBarry Smith    Output Parameter:
800d763cef2SBarry Smith .  v - the vector containing the solution
801d763cef2SBarry Smith 
802d763cef2SBarry Smith    Level: intermediate
803d763cef2SBarry Smith 
804d763cef2SBarry Smith .seealso: TSGetTimeStep()
805d763cef2SBarry Smith 
806d763cef2SBarry Smith .keywords: TS, timestep, get, solution
807d763cef2SBarry Smith @*/
808dfbe8321SBarry Smith PetscErrorCode TSGetSolution(TS ts,Vec *v)
809d763cef2SBarry Smith {
810d763cef2SBarry Smith   PetscFunctionBegin;
8114482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
8124482741eSBarry Smith   PetscValidPointer(v,2);
813d763cef2SBarry Smith   *v = ts->vec_sol_always;
814d763cef2SBarry Smith   PetscFunctionReturn(0);
815d763cef2SBarry Smith }
816d763cef2SBarry Smith 
817bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
8184a2ae208SSatish Balay #undef __FUNCT__
819bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
820d763cef2SBarry Smith /*@C
821bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
822d763cef2SBarry Smith 
823bdad233fSMatthew Knepley   Not collective
824d763cef2SBarry Smith 
825bdad233fSMatthew Knepley   Input Parameters:
826bdad233fSMatthew Knepley + ts   - The TS
827bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
828d763cef2SBarry Smith .vb
829d763cef2SBarry Smith          U_t = A U
830d763cef2SBarry Smith          U_t = A(t) U
831d763cef2SBarry Smith          U_t = F(t,U)
832d763cef2SBarry Smith .ve
833d763cef2SBarry Smith 
834d763cef2SBarry Smith    Level: beginner
835d763cef2SBarry Smith 
836bdad233fSMatthew Knepley .keywords: TS, problem type
837bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
838d763cef2SBarry Smith @*/
839a7cc72afSBarry Smith PetscErrorCode TSSetProblemType(TS ts, TSProblemType type)
840a7cc72afSBarry Smith {
841d763cef2SBarry Smith   PetscFunctionBegin;
8424482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
843bdad233fSMatthew Knepley   ts->problem_type = type;
844d763cef2SBarry Smith   PetscFunctionReturn(0);
845d763cef2SBarry Smith }
846d763cef2SBarry Smith 
847bdad233fSMatthew Knepley #undef __FUNCT__
848bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
849bdad233fSMatthew Knepley /*@C
850bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
851bdad233fSMatthew Knepley 
852bdad233fSMatthew Knepley   Not collective
853bdad233fSMatthew Knepley 
854bdad233fSMatthew Knepley   Input Parameter:
855bdad233fSMatthew Knepley . ts   - The TS
856bdad233fSMatthew Knepley 
857bdad233fSMatthew Knepley   Output Parameter:
858bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
859bdad233fSMatthew Knepley .vb
860bdad233fSMatthew Knepley          U_t = A U
861bdad233fSMatthew Knepley          U_t = A(t) U
862bdad233fSMatthew Knepley          U_t = F(t,U)
863bdad233fSMatthew Knepley .ve
864bdad233fSMatthew Knepley 
865bdad233fSMatthew Knepley    Level: beginner
866bdad233fSMatthew Knepley 
867bdad233fSMatthew Knepley .keywords: TS, problem type
868bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
869bdad233fSMatthew Knepley @*/
870a7cc72afSBarry Smith PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type)
871a7cc72afSBarry Smith {
872bdad233fSMatthew Knepley   PetscFunctionBegin;
8734482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
8744482741eSBarry Smith   PetscValidIntPointer(type,2);
875bdad233fSMatthew Knepley   *type = ts->problem_type;
876bdad233fSMatthew Knepley   PetscFunctionReturn(0);
877bdad233fSMatthew Knepley }
878d763cef2SBarry Smith 
8794a2ae208SSatish Balay #undef __FUNCT__
8804a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
881d763cef2SBarry Smith /*@
882d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
883d763cef2SBarry Smith    of a timestepper.
884d763cef2SBarry Smith 
885d763cef2SBarry Smith    Collective on TS
886d763cef2SBarry Smith 
887d763cef2SBarry Smith    Input Parameter:
888d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
889d763cef2SBarry Smith 
890d763cef2SBarry Smith    Notes:
891d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
892d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
893d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
894d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
895d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
896d763cef2SBarry Smith 
897d763cef2SBarry Smith    Level: advanced
898d763cef2SBarry Smith 
899d763cef2SBarry Smith .keywords: TS, timestep, setup
900d763cef2SBarry Smith 
901d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
902d763cef2SBarry Smith @*/
903dfbe8321SBarry Smith PetscErrorCode TSSetUp(TS ts)
904d763cef2SBarry Smith {
905dfbe8321SBarry Smith   PetscErrorCode ierr;
906d763cef2SBarry Smith 
907d763cef2SBarry Smith   PetscFunctionBegin;
9084482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
90929bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
910d763cef2SBarry Smith   if (!ts->type_name) {
911d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
912d763cef2SBarry Smith   }
913000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
914d763cef2SBarry Smith   ts->setupcalled = 1;
915d763cef2SBarry Smith   PetscFunctionReturn(0);
916d763cef2SBarry Smith }
917d763cef2SBarry Smith 
9184a2ae208SSatish Balay #undef __FUNCT__
9194a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
920d763cef2SBarry Smith /*@C
921d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
922d763cef2SBarry Smith    with TSCreate().
923d763cef2SBarry Smith 
924d763cef2SBarry Smith    Collective on TS
925d763cef2SBarry Smith 
926d763cef2SBarry Smith    Input Parameter:
927d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
928d763cef2SBarry Smith 
929d763cef2SBarry Smith    Level: beginner
930d763cef2SBarry Smith 
931d763cef2SBarry Smith .keywords: TS, timestepper, destroy
932d763cef2SBarry Smith 
933d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
934d763cef2SBarry Smith @*/
935dfbe8321SBarry Smith PetscErrorCode TSDestroy(TS ts)
936d763cef2SBarry Smith {
9376849ba73SBarry Smith   PetscErrorCode ierr;
938a7cc72afSBarry Smith   PetscInt       i;
939d763cef2SBarry Smith 
940d763cef2SBarry Smith   PetscFunctionBegin;
9414482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
942d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
943d763cef2SBarry Smith 
944be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
9450f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
946be0abb6dSBarry Smith 
94794b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
948d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
949000e7ae3SMatthew Knepley   ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);
950329f5518SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
951329f5518SBarry Smith     if (ts->mdestroy[i]) {
952329f5518SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
953329f5518SBarry Smith     }
954329f5518SBarry Smith   }
955*a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
956d763cef2SBarry Smith   PetscFunctionReturn(0);
957d763cef2SBarry Smith }
958d763cef2SBarry Smith 
9594a2ae208SSatish Balay #undef __FUNCT__
9604a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
961d763cef2SBarry Smith /*@C
962d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
963d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
964d763cef2SBarry Smith 
965d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
966d763cef2SBarry Smith 
967d763cef2SBarry Smith    Input Parameter:
968d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
969d763cef2SBarry Smith 
970d763cef2SBarry Smith    Output Parameter:
971d763cef2SBarry Smith .  snes - the nonlinear solver context
972d763cef2SBarry Smith 
973d763cef2SBarry Smith    Notes:
974d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
975d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
97694b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
977d763cef2SBarry Smith 
978d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
979d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
980d763cef2SBarry Smith 
981d763cef2SBarry Smith    Level: beginner
982d763cef2SBarry Smith 
983d763cef2SBarry Smith .keywords: timestep, get, SNES
984d763cef2SBarry Smith @*/
985dfbe8321SBarry Smith PetscErrorCode TSGetSNES(TS ts,SNES *snes)
986d763cef2SBarry Smith {
987d763cef2SBarry Smith   PetscFunctionBegin;
9884482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
9894482741eSBarry Smith   PetscValidPointer(snes,2);
99094b7f48cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
991d763cef2SBarry Smith   *snes = ts->snes;
992d763cef2SBarry Smith   PetscFunctionReturn(0);
993d763cef2SBarry Smith }
994d763cef2SBarry Smith 
9954a2ae208SSatish Balay #undef __FUNCT__
99694b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
997d763cef2SBarry Smith /*@C
99894b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
999d763cef2SBarry Smith    a TS (timestepper) context.
1000d763cef2SBarry Smith 
100194b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1002d763cef2SBarry Smith 
1003d763cef2SBarry Smith    Input Parameter:
1004d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1005d763cef2SBarry Smith 
1006d763cef2SBarry Smith    Output Parameter:
100794b7f48cSBarry Smith .  ksp - the nonlinear solver context
1008d763cef2SBarry Smith 
1009d763cef2SBarry Smith    Notes:
101094b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1011d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1012d763cef2SBarry Smith    KSP and PC contexts as well.
1013d763cef2SBarry Smith 
101494b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
101594b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1016d763cef2SBarry Smith 
1017d763cef2SBarry Smith    Level: beginner
1018d763cef2SBarry Smith 
101994b7f48cSBarry Smith .keywords: timestep, get, KSP
1020d763cef2SBarry Smith @*/
1021dfbe8321SBarry Smith PetscErrorCode TSGetKSP(TS ts,KSP *ksp)
1022d763cef2SBarry Smith {
1023d763cef2SBarry Smith   PetscFunctionBegin;
10244482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
10254482741eSBarry Smith   PetscValidPointer(ksp,2);
102629bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
102794b7f48cSBarry Smith   *ksp = ts->ksp;
1028d763cef2SBarry Smith   PetscFunctionReturn(0);
1029d763cef2SBarry Smith }
1030d763cef2SBarry Smith 
1031d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1032d763cef2SBarry Smith 
10334a2ae208SSatish Balay #undef __FUNCT__
1034adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1035adb62b0dSMatthew Knepley /*@
1036adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1037adb62b0dSMatthew Knepley    maximum time for iteration.
1038adb62b0dSMatthew Knepley 
1039adb62b0dSMatthew Knepley    Collective on TS
1040adb62b0dSMatthew Knepley 
1041adb62b0dSMatthew Knepley    Input Parameters:
1042adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1043adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1044adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1045adb62b0dSMatthew Knepley 
1046adb62b0dSMatthew Knepley    Level: intermediate
1047adb62b0dSMatthew Knepley 
1048adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1049adb62b0dSMatthew Knepley @*/
1050a7cc72afSBarry Smith PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1051adb62b0dSMatthew Knepley {
1052adb62b0dSMatthew Knepley   PetscFunctionBegin;
10534482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1054abc0a331SBarry Smith   if (maxsteps) {
10554482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1056adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1057adb62b0dSMatthew Knepley   }
1058abc0a331SBarry Smith   if (maxtime ) {
10594482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1060adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1061adb62b0dSMatthew Knepley   }
1062adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1063adb62b0dSMatthew Knepley }
1064adb62b0dSMatthew Knepley 
1065adb62b0dSMatthew Knepley #undef __FUNCT__
10664a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1067d763cef2SBarry Smith /*@
1068d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1069d763cef2SBarry Smith    maximum time for iteration.
1070d763cef2SBarry Smith 
1071d763cef2SBarry Smith    Collective on TS
1072d763cef2SBarry Smith 
1073d763cef2SBarry Smith    Input Parameters:
1074d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1075d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1076d763cef2SBarry Smith -  maxtime - final time to iterate to
1077d763cef2SBarry Smith 
1078d763cef2SBarry Smith    Options Database Keys:
1079d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1080d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1081d763cef2SBarry Smith 
1082d763cef2SBarry Smith    Notes:
1083d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1084d763cef2SBarry Smith 
1085d763cef2SBarry Smith    Level: intermediate
1086d763cef2SBarry Smith 
1087d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1088d763cef2SBarry Smith @*/
1089a7cc72afSBarry Smith PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1090d763cef2SBarry Smith {
1091d763cef2SBarry Smith   PetscFunctionBegin;
10924482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1093d763cef2SBarry Smith   ts->max_steps = maxsteps;
1094d763cef2SBarry Smith   ts->max_time  = maxtime;
1095d763cef2SBarry Smith   PetscFunctionReturn(0);
1096d763cef2SBarry Smith }
1097d763cef2SBarry Smith 
10984a2ae208SSatish Balay #undef __FUNCT__
10994a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1100d763cef2SBarry Smith /*@
1101d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1102d763cef2SBarry Smith    for use by the TS routines.
1103d763cef2SBarry Smith 
1104d763cef2SBarry Smith    Collective on TS and Vec
1105d763cef2SBarry Smith 
1106d763cef2SBarry Smith    Input Parameters:
1107d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1108d763cef2SBarry Smith -  x - the solution vector
1109d763cef2SBarry Smith 
1110d763cef2SBarry Smith    Level: beginner
1111d763cef2SBarry Smith 
1112d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1113d763cef2SBarry Smith @*/
1114dfbe8321SBarry Smith PetscErrorCode TSSetSolution(TS ts,Vec x)
1115d763cef2SBarry Smith {
1116d763cef2SBarry Smith   PetscFunctionBegin;
11174482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
11184482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1119d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1120d763cef2SBarry Smith   PetscFunctionReturn(0);
1121d763cef2SBarry Smith }
1122d763cef2SBarry Smith 
1123e74ef692SMatthew Knepley #undef __FUNCT__
1124e74ef692SMatthew Knepley #define __FUNCT__ "TSSetRhsBC"
1125ac226902SBarry Smith /*@C
1126000e7ae3SMatthew Knepley   TSSetRhsBC - Sets the function which applies boundary conditions
1127000e7ae3SMatthew Knepley   to the Rhs of each system.
1128000e7ae3SMatthew Knepley 
1129000e7ae3SMatthew Knepley   Collective on TS
1130000e7ae3SMatthew Knepley 
1131000e7ae3SMatthew Knepley   Input Parameters:
1132000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1133000e7ae3SMatthew Knepley - func - The function
1134000e7ae3SMatthew Knepley 
1135000e7ae3SMatthew Knepley   Calling sequence of func:
1136000e7ae3SMatthew Knepley . func (TS ts, Vec rhs, void *ctx);
1137000e7ae3SMatthew Knepley 
1138000e7ae3SMatthew Knepley + rhs - The current rhs vector
1139000e7ae3SMatthew Knepley - ctx - The user-context
1140000e7ae3SMatthew Knepley 
1141000e7ae3SMatthew Knepley   Level: intermediate
1142000e7ae3SMatthew Knepley 
1143000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1144000e7ae3SMatthew Knepley @*/
11456849ba73SBarry Smith PetscErrorCode TSSetRhsBC(TS ts, PetscErrorCode (*func)(TS, Vec, void *))
1146000e7ae3SMatthew Knepley {
1147000e7ae3SMatthew Knepley   PetscFunctionBegin;
11484482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1149000e7ae3SMatthew Knepley   ts->ops->applyrhsbc = func;
1150000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1151000e7ae3SMatthew Knepley }
1152000e7ae3SMatthew Knepley 
1153e74ef692SMatthew Knepley #undef __FUNCT__
1154e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultRhsBC"
1155000e7ae3SMatthew Knepley /*@
1156000e7ae3SMatthew Knepley   TSDefaultRhsBC - The default boundary condition function which does nothing.
1157000e7ae3SMatthew Knepley 
1158000e7ae3SMatthew Knepley   Collective on TS
1159000e7ae3SMatthew Knepley 
1160000e7ae3SMatthew Knepley   Input Parameters:
1161000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1162000e7ae3SMatthew Knepley . rhs - The Rhs
1163000e7ae3SMatthew Knepley - ctx - The user-context
1164000e7ae3SMatthew Knepley 
1165000e7ae3SMatthew Knepley   Level: developer
1166000e7ae3SMatthew Knepley 
1167000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1168000e7ae3SMatthew Knepley @*/
1169dfbe8321SBarry Smith PetscErrorCode TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1170000e7ae3SMatthew Knepley {
1171000e7ae3SMatthew Knepley   PetscFunctionBegin;
1172000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1173000e7ae3SMatthew Knepley }
1174000e7ae3SMatthew Knepley 
1175e74ef692SMatthew Knepley #undef __FUNCT__
1176e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSystemMatrixBC"
1177ac226902SBarry Smith /*@C
1178000e7ae3SMatthew Knepley   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1179000e7ae3SMatthew Knepley   to the system matrix and preconditioner of each system.
1180000e7ae3SMatthew Knepley 
1181000e7ae3SMatthew Knepley   Collective on TS
1182000e7ae3SMatthew Knepley 
1183000e7ae3SMatthew Knepley   Input Parameters:
1184000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1185000e7ae3SMatthew Knepley - func - The function
1186000e7ae3SMatthew Knepley 
1187000e7ae3SMatthew Knepley   Calling sequence of func:
1188000e7ae3SMatthew Knepley . func (TS ts, Mat A, Mat B, void *ctx);
1189000e7ae3SMatthew Knepley 
1190000e7ae3SMatthew Knepley + A   - The current system matrix
1191000e7ae3SMatthew Knepley . B   - The current preconditioner
1192000e7ae3SMatthew Knepley - ctx - The user-context
1193000e7ae3SMatthew Knepley 
1194000e7ae3SMatthew Knepley   Level: intermediate
1195000e7ae3SMatthew Knepley 
1196000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1197000e7ae3SMatthew Knepley @*/
11986849ba73SBarry Smith PetscErrorCode TSSetSystemMatrixBC(TS ts, PetscErrorCode (*func)(TS, Mat, Mat, void *))
1199000e7ae3SMatthew Knepley {
1200000e7ae3SMatthew Knepley   PetscFunctionBegin;
12014482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1202000e7ae3SMatthew Knepley   ts->ops->applymatrixbc = func;
1203000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1204000e7ae3SMatthew Knepley }
1205000e7ae3SMatthew Knepley 
1206e74ef692SMatthew Knepley #undef __FUNCT__
1207e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSystemMatrixBC"
1208000e7ae3SMatthew Knepley /*@
1209000e7ae3SMatthew Knepley   TSDefaultSystemMatrixBC - The default boundary condition function which
1210000e7ae3SMatthew Knepley   does nothing.
1211000e7ae3SMatthew Knepley 
1212000e7ae3SMatthew Knepley   Collective on TS
1213000e7ae3SMatthew Knepley 
1214000e7ae3SMatthew Knepley   Input Parameters:
1215000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1216000e7ae3SMatthew Knepley . A   - The system matrix
1217000e7ae3SMatthew Knepley . B   - The preconditioner
1218000e7ae3SMatthew Knepley - ctx - The user-context
1219000e7ae3SMatthew Knepley 
1220000e7ae3SMatthew Knepley   Level: developer
1221000e7ae3SMatthew Knepley 
1222000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1223000e7ae3SMatthew Knepley @*/
1224dfbe8321SBarry Smith PetscErrorCode TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1225000e7ae3SMatthew Knepley {
1226000e7ae3SMatthew Knepley   PetscFunctionBegin;
1227000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1228000e7ae3SMatthew Knepley }
1229000e7ae3SMatthew Knepley 
1230e74ef692SMatthew Knepley #undef __FUNCT__
1231e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSolutionBC"
1232ac226902SBarry Smith /*@C
1233000e7ae3SMatthew Knepley   TSSetSolutionBC - Sets the function which applies boundary conditions
1234000e7ae3SMatthew Knepley   to the solution of each system. This is necessary in nonlinear systems
1235000e7ae3SMatthew Knepley   which time dependent boundary conditions.
1236000e7ae3SMatthew Knepley 
1237000e7ae3SMatthew Knepley   Collective on TS
1238000e7ae3SMatthew Knepley 
1239000e7ae3SMatthew Knepley   Input Parameters:
1240000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1241000e7ae3SMatthew Knepley - func - The function
1242000e7ae3SMatthew Knepley 
1243000e7ae3SMatthew Knepley   Calling sequence of func:
1244000e7ae3SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
1245000e7ae3SMatthew Knepley 
1246000e7ae3SMatthew Knepley + sol - The current solution vector
1247000e7ae3SMatthew Knepley - ctx - The user-context
1248000e7ae3SMatthew Knepley 
1249000e7ae3SMatthew Knepley   Level: intermediate
1250000e7ae3SMatthew Knepley 
1251000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1252000e7ae3SMatthew Knepley @*/
12536849ba73SBarry Smith PetscErrorCode TSSetSolutionBC(TS ts, PetscErrorCode (*func)(TS, Vec, void *))
1254000e7ae3SMatthew Knepley {
1255000e7ae3SMatthew Knepley   PetscFunctionBegin;
12564482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1257000e7ae3SMatthew Knepley   ts->ops->applysolbc = func;
1258000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1259000e7ae3SMatthew Knepley }
1260000e7ae3SMatthew Knepley 
1261e74ef692SMatthew Knepley #undef __FUNCT__
1262e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSolutionBC"
1263000e7ae3SMatthew Knepley /*@
1264000e7ae3SMatthew Knepley   TSDefaultSolutionBC - The default boundary condition function which
1265000e7ae3SMatthew Knepley   does nothing.
1266000e7ae3SMatthew Knepley 
1267000e7ae3SMatthew Knepley   Collective on TS
1268000e7ae3SMatthew Knepley 
1269000e7ae3SMatthew Knepley   Input Parameters:
1270000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1271000e7ae3SMatthew Knepley . sol - The solution
1272000e7ae3SMatthew Knepley - ctx - The user-context
1273000e7ae3SMatthew Knepley 
1274000e7ae3SMatthew Knepley   Level: developer
1275000e7ae3SMatthew Knepley 
1276000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1277000e7ae3SMatthew Knepley @*/
1278dfbe8321SBarry Smith PetscErrorCode TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1279000e7ae3SMatthew Knepley {
1280000e7ae3SMatthew Knepley   PetscFunctionBegin;
1281000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1282000e7ae3SMatthew Knepley }
1283000e7ae3SMatthew Knepley 
1284e74ef692SMatthew Knepley #undef __FUNCT__
1285e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1286ac226902SBarry Smith /*@C
1287000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1288000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1289000e7ae3SMatthew Knepley 
1290000e7ae3SMatthew Knepley   Collective on TS
1291000e7ae3SMatthew Knepley 
1292000e7ae3SMatthew Knepley   Input Parameters:
1293000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1294000e7ae3SMatthew Knepley - func - The function
1295000e7ae3SMatthew Knepley 
1296000e7ae3SMatthew Knepley   Calling sequence of func:
1297000e7ae3SMatthew Knepley . func (TS ts);
1298000e7ae3SMatthew Knepley 
1299000e7ae3SMatthew Knepley   Level: intermediate
1300000e7ae3SMatthew Knepley 
1301000e7ae3SMatthew Knepley .keywords: TS, timestep
1302000e7ae3SMatthew Knepley @*/
13036849ba73SBarry Smith PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1304000e7ae3SMatthew Knepley {
1305000e7ae3SMatthew Knepley   PetscFunctionBegin;
13064482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1307000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1308000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1309000e7ae3SMatthew Knepley }
1310000e7ae3SMatthew Knepley 
1311e74ef692SMatthew Knepley #undef __FUNCT__
1312e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1313000e7ae3SMatthew Knepley /*@
1314000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1315000e7ae3SMatthew Knepley 
1316000e7ae3SMatthew Knepley   Collective on TS
1317000e7ae3SMatthew Knepley 
1318000e7ae3SMatthew Knepley   Input Parameters:
1319000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1320000e7ae3SMatthew Knepley 
1321000e7ae3SMatthew Knepley   Level: developer
1322000e7ae3SMatthew Knepley 
1323000e7ae3SMatthew Knepley .keywords: TS, timestep
1324000e7ae3SMatthew Knepley @*/
1325dfbe8321SBarry Smith PetscErrorCode TSDefaultPreStep(TS ts)
1326000e7ae3SMatthew Knepley {
1327000e7ae3SMatthew Knepley   PetscFunctionBegin;
1328000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1329000e7ae3SMatthew Knepley }
1330000e7ae3SMatthew Knepley 
1331e74ef692SMatthew Knepley #undef __FUNCT__
1332e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1333ac226902SBarry Smith /*@C
1334000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1335000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1336000e7ae3SMatthew Knepley   the time step.
1337000e7ae3SMatthew Knepley 
1338000e7ae3SMatthew Knepley   Collective on TS
1339000e7ae3SMatthew Knepley 
1340000e7ae3SMatthew Knepley   Input Parameters:
1341000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1342000e7ae3SMatthew Knepley - func - The function
1343000e7ae3SMatthew Knepley 
1344000e7ae3SMatthew Knepley   Calling sequence of func:
1345000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1346000e7ae3SMatthew Knepley 
1347000e7ae3SMatthew Knepley + t   - The current time
1348000e7ae3SMatthew Knepley - dt  - The current time step
1349000e7ae3SMatthew Knepley 
1350000e7ae3SMatthew Knepley   Level: intermediate
1351000e7ae3SMatthew Knepley 
1352000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1353000e7ae3SMatthew Knepley @*/
13546849ba73SBarry Smith PetscErrorCode TSSetUpdate(TS ts, PetscErrorCode (*func)(TS, PetscReal, PetscReal *))
1355000e7ae3SMatthew Knepley {
1356000e7ae3SMatthew Knepley   PetscFunctionBegin;
13574482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1358000e7ae3SMatthew Knepley   ts->ops->update = func;
1359000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1360000e7ae3SMatthew Knepley }
1361000e7ae3SMatthew Knepley 
1362e74ef692SMatthew Knepley #undef __FUNCT__
1363e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1364000e7ae3SMatthew Knepley /*@
1365000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1366000e7ae3SMatthew Knepley 
1367000e7ae3SMatthew Knepley   Collective on TS
1368000e7ae3SMatthew Knepley 
1369000e7ae3SMatthew Knepley   Input Parameters:
1370000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1371000e7ae3SMatthew Knepley - t   - The current time
1372000e7ae3SMatthew Knepley 
1373000e7ae3SMatthew Knepley   Output Parameters:
1374000e7ae3SMatthew Knepley . dt  - The current time step
1375000e7ae3SMatthew Knepley 
1376000e7ae3SMatthew Knepley   Level: developer
1377000e7ae3SMatthew Knepley 
1378000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1379000e7ae3SMatthew Knepley @*/
1380dfbe8321SBarry Smith PetscErrorCode TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1381000e7ae3SMatthew Knepley {
1382000e7ae3SMatthew Knepley   PetscFunctionBegin;
1383000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1384000e7ae3SMatthew Knepley }
1385000e7ae3SMatthew Knepley 
1386e74ef692SMatthew Knepley #undef __FUNCT__
1387e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1388ac226902SBarry Smith /*@C
1389000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1390000e7ae3SMatthew Knepley   called once at the end of time stepping.
1391000e7ae3SMatthew Knepley 
1392000e7ae3SMatthew Knepley   Collective on TS
1393000e7ae3SMatthew Knepley 
1394000e7ae3SMatthew Knepley   Input Parameters:
1395000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1396000e7ae3SMatthew Knepley - func - The function
1397000e7ae3SMatthew Knepley 
1398000e7ae3SMatthew Knepley   Calling sequence of func:
1399000e7ae3SMatthew Knepley . func (TS ts);
1400000e7ae3SMatthew Knepley 
1401000e7ae3SMatthew Knepley   Level: intermediate
1402000e7ae3SMatthew Knepley 
1403000e7ae3SMatthew Knepley .keywords: TS, timestep
1404000e7ae3SMatthew Knepley @*/
14056849ba73SBarry Smith PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1406000e7ae3SMatthew Knepley {
1407000e7ae3SMatthew Knepley   PetscFunctionBegin;
14084482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1409000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1410000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1411000e7ae3SMatthew Knepley }
1412000e7ae3SMatthew Knepley 
1413e74ef692SMatthew Knepley #undef __FUNCT__
1414e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1415000e7ae3SMatthew Knepley /*@
1416000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1417000e7ae3SMatthew Knepley 
1418000e7ae3SMatthew Knepley   Collective on TS
1419000e7ae3SMatthew Knepley 
1420000e7ae3SMatthew Knepley   Input Parameters:
1421000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1422000e7ae3SMatthew Knepley 
1423000e7ae3SMatthew Knepley   Level: developer
1424000e7ae3SMatthew Knepley 
1425000e7ae3SMatthew Knepley .keywords: TS, timestep
1426000e7ae3SMatthew Knepley @*/
1427dfbe8321SBarry Smith PetscErrorCode TSDefaultPostStep(TS ts)
1428000e7ae3SMatthew Knepley {
1429000e7ae3SMatthew Knepley   PetscFunctionBegin;
1430000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1431000e7ae3SMatthew Knepley }
1432000e7ae3SMatthew Knepley 
1433d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1434d763cef2SBarry Smith 
14354a2ae208SSatish Balay #undef __FUNCT__
14364a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
1437d763cef2SBarry Smith /*@C
1438d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1439d763cef2SBarry Smith    timestep to display the iteration's  progress.
1440d763cef2SBarry Smith 
1441d763cef2SBarry Smith    Collective on TS
1442d763cef2SBarry Smith 
1443d763cef2SBarry Smith    Input Parameters:
1444d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1445d763cef2SBarry Smith .  func - monitoring routine
1446329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1447b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1448b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1449b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1450d763cef2SBarry Smith 
1451d763cef2SBarry Smith    Calling sequence of func:
1452a7cc72afSBarry Smith $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1453d763cef2SBarry Smith 
1454d763cef2SBarry Smith +    ts - the TS context
1455d763cef2SBarry Smith .    steps - iteration number
14561f06c33eSBarry Smith .    time - current time
1457d763cef2SBarry Smith .    x - current iterate
1458d763cef2SBarry Smith -    mctx - [optional] monitoring context
1459d763cef2SBarry Smith 
1460d763cef2SBarry Smith    Notes:
1461d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1462d763cef2SBarry Smith    already has been loaded.
1463d763cef2SBarry Smith 
1464d763cef2SBarry Smith    Level: intermediate
1465d763cef2SBarry Smith 
1466d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1467d763cef2SBarry Smith 
1468d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
1469d763cef2SBarry Smith @*/
1470a7cc72afSBarry Smith PetscErrorCode TSSetMonitor(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1471d763cef2SBarry Smith {
1472d763cef2SBarry Smith   PetscFunctionBegin;
14734482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1474d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
147529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1476d763cef2SBarry Smith   }
1477d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1478329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1479d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1480d763cef2SBarry Smith   PetscFunctionReturn(0);
1481d763cef2SBarry Smith }
1482d763cef2SBarry Smith 
14834a2ae208SSatish Balay #undef __FUNCT__
14844a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
1485d763cef2SBarry Smith /*@C
1486d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1487d763cef2SBarry Smith 
1488d763cef2SBarry Smith    Collective on TS
1489d763cef2SBarry Smith 
1490d763cef2SBarry Smith    Input Parameters:
1491d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1492d763cef2SBarry Smith 
1493d763cef2SBarry Smith    Notes:
1494d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1495d763cef2SBarry Smith 
1496d763cef2SBarry Smith    Level: intermediate
1497d763cef2SBarry Smith 
1498d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1499d763cef2SBarry Smith 
1500d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
1501d763cef2SBarry Smith @*/
1502dfbe8321SBarry Smith PetscErrorCode TSClearMonitor(TS ts)
1503d763cef2SBarry Smith {
1504d763cef2SBarry Smith   PetscFunctionBegin;
15054482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1506d763cef2SBarry Smith   ts->numbermonitors = 0;
1507d763cef2SBarry Smith   PetscFunctionReturn(0);
1508d763cef2SBarry Smith }
1509d763cef2SBarry Smith 
15104a2ae208SSatish Balay #undef __FUNCT__
15114a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
1512a7cc72afSBarry Smith PetscErrorCode TSDefaultMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1513d763cef2SBarry Smith {
1514dfbe8321SBarry Smith   PetscErrorCode ierr;
1515d132466eSBarry Smith 
1516d763cef2SBarry Smith   PetscFunctionBegin;
151777431f27SBarry Smith   ierr = PetscPrintf(ts->comm,"timestep %D dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1518d763cef2SBarry Smith   PetscFunctionReturn(0);
1519d763cef2SBarry Smith }
1520d763cef2SBarry Smith 
15214a2ae208SSatish Balay #undef __FUNCT__
15224a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1523d763cef2SBarry Smith /*@
1524d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1525d763cef2SBarry Smith 
1526d763cef2SBarry Smith    Collective on TS
1527d763cef2SBarry Smith 
1528d763cef2SBarry Smith    Input Parameter:
1529d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1530d763cef2SBarry Smith 
1531d763cef2SBarry Smith    Output Parameters:
1532d763cef2SBarry Smith +  steps - number of iterations until termination
1533142b95e3SSatish Balay -  ptime - time until termination
1534d763cef2SBarry Smith 
1535d763cef2SBarry Smith    Level: beginner
1536d763cef2SBarry Smith 
1537d763cef2SBarry Smith .keywords: TS, timestep, solve
1538d763cef2SBarry Smith 
1539d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1540d763cef2SBarry Smith @*/
1541a7cc72afSBarry Smith PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1542d763cef2SBarry Smith {
1543dfbe8321SBarry Smith   PetscErrorCode ierr;
1544d763cef2SBarry Smith 
1545d763cef2SBarry Smith   PetscFunctionBegin;
15464482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1547d405a339SMatthew Knepley   if (!ts->setupcalled) {
1548d405a339SMatthew Knepley     ierr = TSSetUp(ts);CHKERRQ(ierr);
1549d405a339SMatthew Knepley   }
1550d405a339SMatthew Knepley 
1551d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1552d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1553000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1554d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1555d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1556d405a339SMatthew Knepley 
15574bb05414SBarry Smith   if (!PetscPreLoadingOn) {
15584bb05414SBarry Smith     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1559d405a339SMatthew Knepley   }
1560d763cef2SBarry Smith   PetscFunctionReturn(0);
1561d763cef2SBarry Smith }
1562d763cef2SBarry Smith 
15634a2ae208SSatish Balay #undef __FUNCT__
15644a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1565d763cef2SBarry Smith /*
1566d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1567d763cef2SBarry Smith */
1568a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1569d763cef2SBarry Smith {
15706849ba73SBarry Smith   PetscErrorCode ierr;
1571a7cc72afSBarry Smith   PetscInt i,n = ts->numbermonitors;
1572d763cef2SBarry Smith 
1573d763cef2SBarry Smith   PetscFunctionBegin;
1574d763cef2SBarry Smith   for (i=0; i<n; i++) {
1575142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1576d763cef2SBarry Smith   }
1577d763cef2SBarry Smith   PetscFunctionReturn(0);
1578d763cef2SBarry Smith }
1579d763cef2SBarry Smith 
1580d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1581d763cef2SBarry Smith 
15824a2ae208SSatish Balay #undef __FUNCT__
15834a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1584d763cef2SBarry Smith /*@C
1585d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1586d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1587d763cef2SBarry Smith 
1588d763cef2SBarry Smith    Collective on TS
1589d763cef2SBarry Smith 
1590d763cef2SBarry Smith    Input Parameters:
1591d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1592d763cef2SBarry Smith .  label - the title to put in the title bar
15937c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1594d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1595d763cef2SBarry Smith 
1596d763cef2SBarry Smith    Output Parameter:
1597d763cef2SBarry Smith .  draw - the drawing context
1598d763cef2SBarry Smith 
1599d763cef2SBarry Smith    Options Database Key:
1600d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1601d763cef2SBarry Smith 
1602d763cef2SBarry Smith    Notes:
1603b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1604d763cef2SBarry Smith 
1605d763cef2SBarry Smith    Level: intermediate
1606d763cef2SBarry Smith 
16077c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1608d763cef2SBarry Smith 
1609d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
16107c922b88SBarry Smith 
1611d763cef2SBarry Smith @*/
1612dfbe8321SBarry Smith PetscErrorCode TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1613d763cef2SBarry Smith {
1614b0a32e0cSBarry Smith   PetscDraw      win;
1615dfbe8321SBarry Smith   PetscErrorCode ierr;
1616d763cef2SBarry Smith 
1617d763cef2SBarry Smith   PetscFunctionBegin;
1618b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1619b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1620b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1621b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1622d763cef2SBarry Smith 
162352e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1624d763cef2SBarry Smith   PetscFunctionReturn(0);
1625d763cef2SBarry Smith }
1626d763cef2SBarry Smith 
16274a2ae208SSatish Balay #undef __FUNCT__
16284a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
1629a7cc72afSBarry Smith PetscErrorCode TSLGMonitor(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1630d763cef2SBarry Smith {
1631b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
163287828ca2SBarry Smith   PetscReal      x,y = ptime;
1633dfbe8321SBarry Smith   PetscErrorCode ierr;
1634d763cef2SBarry Smith 
1635d763cef2SBarry Smith   PetscFunctionBegin;
16367c922b88SBarry Smith   if (!monctx) {
16377c922b88SBarry Smith     MPI_Comm    comm;
1638b0a32e0cSBarry Smith     PetscViewer viewer;
16397c922b88SBarry Smith 
16407c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1641b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1642b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
16437c922b88SBarry Smith   }
16447c922b88SBarry Smith 
1645b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
164687828ca2SBarry Smith   x = (PetscReal)n;
1647b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1648d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1649b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1650d763cef2SBarry Smith   }
1651d763cef2SBarry Smith   PetscFunctionReturn(0);
1652d763cef2SBarry Smith }
1653d763cef2SBarry Smith 
16544a2ae208SSatish Balay #undef __FUNCT__
16554a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1656d763cef2SBarry Smith /*@C
1657d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1658d763cef2SBarry Smith    with TSLGMonitorCreate().
1659d763cef2SBarry Smith 
1660b0a32e0cSBarry Smith    Collective on PetscDrawLG
1661d763cef2SBarry Smith 
1662d763cef2SBarry Smith    Input Parameter:
1663d763cef2SBarry Smith .  draw - the drawing context
1664d763cef2SBarry Smith 
1665d763cef2SBarry Smith    Level: intermediate
1666d763cef2SBarry Smith 
1667d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1668d763cef2SBarry Smith 
1669d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1670d763cef2SBarry Smith @*/
1671dfbe8321SBarry Smith PetscErrorCode TSLGMonitorDestroy(PetscDrawLG drawlg)
1672d763cef2SBarry Smith {
1673b0a32e0cSBarry Smith   PetscDraw      draw;
1674dfbe8321SBarry Smith   PetscErrorCode ierr;
1675d763cef2SBarry Smith 
1676d763cef2SBarry Smith   PetscFunctionBegin;
1677b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1678b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1679b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1680d763cef2SBarry Smith   PetscFunctionReturn(0);
1681d763cef2SBarry Smith }
1682d763cef2SBarry Smith 
16834a2ae208SSatish Balay #undef __FUNCT__
16844a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1685d763cef2SBarry Smith /*@
1686d763cef2SBarry Smith    TSGetTime - Gets the current time.
1687d763cef2SBarry Smith 
1688d763cef2SBarry Smith    Not Collective
1689d763cef2SBarry Smith 
1690d763cef2SBarry Smith    Input Parameter:
1691d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1692d763cef2SBarry Smith 
1693d763cef2SBarry Smith    Output Parameter:
1694d763cef2SBarry Smith .  t  - the current time
1695d763cef2SBarry Smith 
1696d763cef2SBarry Smith    Contributed by: Matthew Knepley
1697d763cef2SBarry Smith 
1698d763cef2SBarry Smith    Level: beginner
1699d763cef2SBarry Smith 
1700d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1701d763cef2SBarry Smith 
1702d763cef2SBarry Smith .keywords: TS, get, time
1703d763cef2SBarry Smith @*/
1704dfbe8321SBarry Smith PetscErrorCode TSGetTime(TS ts,PetscReal* t)
1705d763cef2SBarry Smith {
1706d763cef2SBarry Smith   PetscFunctionBegin;
17074482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
17084482741eSBarry Smith   PetscValidDoublePointer(t,2);
1709d763cef2SBarry Smith   *t = ts->ptime;
1710d763cef2SBarry Smith   PetscFunctionReturn(0);
1711d763cef2SBarry Smith }
1712d763cef2SBarry Smith 
17134a2ae208SSatish Balay #undef __FUNCT__
17144a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1715d763cef2SBarry Smith /*@C
1716d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1717d763cef2SBarry Smith    TS options in the database.
1718d763cef2SBarry Smith 
1719d763cef2SBarry Smith    Collective on TS
1720d763cef2SBarry Smith 
1721d763cef2SBarry Smith    Input Parameter:
1722d763cef2SBarry Smith +  ts     - The TS context
1723d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1724d763cef2SBarry Smith 
1725d763cef2SBarry Smith    Notes:
1726d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1727d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1728d763cef2SBarry Smith    hyphen.
1729d763cef2SBarry Smith 
1730d763cef2SBarry Smith    Contributed by: Matthew Knepley
1731d763cef2SBarry Smith 
1732d763cef2SBarry Smith    Level: advanced
1733d763cef2SBarry Smith 
1734d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1735d763cef2SBarry Smith 
1736d763cef2SBarry Smith .seealso: TSSetFromOptions()
1737d763cef2SBarry Smith 
1738d763cef2SBarry Smith @*/
1739dfbe8321SBarry Smith PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[])
1740d763cef2SBarry Smith {
1741dfbe8321SBarry Smith   PetscErrorCode ierr;
1742d763cef2SBarry Smith 
1743d763cef2SBarry Smith   PetscFunctionBegin;
17444482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1745d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1746d763cef2SBarry Smith   switch(ts->problem_type) {
1747d763cef2SBarry Smith     case TS_NONLINEAR:
1748d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1749d763cef2SBarry Smith       break;
1750d763cef2SBarry Smith     case TS_LINEAR:
175194b7f48cSBarry Smith       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1752d763cef2SBarry Smith       break;
1753d763cef2SBarry Smith   }
1754d763cef2SBarry Smith   PetscFunctionReturn(0);
1755d763cef2SBarry Smith }
1756d763cef2SBarry Smith 
1757d763cef2SBarry Smith 
17584a2ae208SSatish Balay #undef __FUNCT__
17594a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1760d763cef2SBarry Smith /*@C
1761d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1762d763cef2SBarry Smith    TS options in the database.
1763d763cef2SBarry Smith 
1764d763cef2SBarry Smith    Collective on TS
1765d763cef2SBarry Smith 
1766d763cef2SBarry Smith    Input Parameter:
1767d763cef2SBarry Smith +  ts     - The TS context
1768d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1769d763cef2SBarry Smith 
1770d763cef2SBarry Smith    Notes:
1771d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1772d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1773d763cef2SBarry Smith    hyphen.
1774d763cef2SBarry Smith 
1775d763cef2SBarry Smith    Contributed by: Matthew Knepley
1776d763cef2SBarry Smith 
1777d763cef2SBarry Smith    Level: advanced
1778d763cef2SBarry Smith 
1779d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1780d763cef2SBarry Smith 
1781d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1782d763cef2SBarry Smith 
1783d763cef2SBarry Smith @*/
1784dfbe8321SBarry Smith PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[])
1785d763cef2SBarry Smith {
1786dfbe8321SBarry Smith   PetscErrorCode ierr;
1787d763cef2SBarry Smith 
1788d763cef2SBarry Smith   PetscFunctionBegin;
17894482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1790d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1791d763cef2SBarry Smith   switch(ts->problem_type) {
1792d763cef2SBarry Smith     case TS_NONLINEAR:
1793d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1794d763cef2SBarry Smith       break;
1795d763cef2SBarry Smith     case TS_LINEAR:
179694b7f48cSBarry Smith       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1797d763cef2SBarry Smith       break;
1798d763cef2SBarry Smith   }
1799d763cef2SBarry Smith   PetscFunctionReturn(0);
1800d763cef2SBarry Smith }
1801d763cef2SBarry Smith 
18024a2ae208SSatish Balay #undef __FUNCT__
18034a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1804d763cef2SBarry Smith /*@C
1805d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1806d763cef2SBarry Smith    TS options in the database.
1807d763cef2SBarry Smith 
1808d763cef2SBarry Smith    Not Collective
1809d763cef2SBarry Smith 
1810d763cef2SBarry Smith    Input Parameter:
1811d763cef2SBarry Smith .  ts - The TS context
1812d763cef2SBarry Smith 
1813d763cef2SBarry Smith    Output Parameter:
1814d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1815d763cef2SBarry Smith 
1816d763cef2SBarry Smith    Contributed by: Matthew Knepley
1817d763cef2SBarry Smith 
1818d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1819d763cef2SBarry Smith    sufficient length to hold the prefix.
1820d763cef2SBarry Smith 
1821d763cef2SBarry Smith    Level: intermediate
1822d763cef2SBarry Smith 
1823d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1824d763cef2SBarry Smith 
1825d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1826d763cef2SBarry Smith @*/
1827dfbe8321SBarry Smith PetscErrorCode TSGetOptionsPrefix(TS ts,char *prefix[])
1828d763cef2SBarry Smith {
1829dfbe8321SBarry Smith   PetscErrorCode ierr;
1830d763cef2SBarry Smith 
1831d763cef2SBarry Smith   PetscFunctionBegin;
18324482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
18334482741eSBarry Smith   PetscValidPointer(prefix,2);
1834d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1835d763cef2SBarry Smith   PetscFunctionReturn(0);
1836d763cef2SBarry Smith }
1837d763cef2SBarry Smith 
18384a2ae208SSatish Balay #undef __FUNCT__
18394a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1840d763cef2SBarry Smith /*@C
1841d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1842d763cef2SBarry Smith 
1843d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1844d763cef2SBarry Smith 
1845d763cef2SBarry Smith    Input Parameter:
1846d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1847d763cef2SBarry Smith 
1848d763cef2SBarry Smith    Output Parameters:
1849d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1850d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1851d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1852d763cef2SBarry Smith 
1853d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1854d763cef2SBarry Smith 
1855d763cef2SBarry Smith    Contributed by: Matthew Knepley
1856d763cef2SBarry Smith 
1857d763cef2SBarry Smith    Level: intermediate
1858d763cef2SBarry Smith 
1859d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1860d763cef2SBarry Smith 
1861d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1862d763cef2SBarry Smith 
1863d763cef2SBarry Smith @*/
1864dfbe8321SBarry Smith PetscErrorCode TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1865d763cef2SBarry Smith {
1866d763cef2SBarry Smith   PetscFunctionBegin;
18674482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1868d763cef2SBarry Smith   if (A)   *A = ts->A;
1869d763cef2SBarry Smith   if (M)   *M = ts->B;
1870d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1871d763cef2SBarry Smith   PetscFunctionReturn(0);
1872d763cef2SBarry Smith }
1873d763cef2SBarry Smith 
18744a2ae208SSatish Balay #undef __FUNCT__
18754a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1876d763cef2SBarry Smith /*@C
1877d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1878d763cef2SBarry Smith 
1879d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1880d763cef2SBarry Smith 
1881d763cef2SBarry Smith    Input Parameter:
1882d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1883d763cef2SBarry Smith 
1884d763cef2SBarry Smith    Output Parameters:
1885d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1886d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1887d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1888d763cef2SBarry Smith 
1889d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1890d763cef2SBarry Smith 
1891d763cef2SBarry Smith    Contributed by: Matthew Knepley
1892d763cef2SBarry Smith 
1893d763cef2SBarry Smith    Level: intermediate
1894d763cef2SBarry Smith 
1895d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1896d763cef2SBarry Smith 
1897d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1898d763cef2SBarry Smith @*/
1899dfbe8321SBarry Smith PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1900d763cef2SBarry Smith {
1901dfbe8321SBarry Smith   PetscErrorCode ierr;
1902d763cef2SBarry Smith 
1903d763cef2SBarry Smith   PetscFunctionBegin;
1904d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1905d763cef2SBarry Smith   PetscFunctionReturn(0);
1906d763cef2SBarry Smith }
1907d763cef2SBarry Smith 
19081713a123SBarry Smith #undef __FUNCT__
19091713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
19101713a123SBarry Smith /*@C
19111713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
19121713a123SBarry Smith    VecView() for the solution at each timestep
19131713a123SBarry Smith 
19141713a123SBarry Smith    Collective on TS
19151713a123SBarry Smith 
19161713a123SBarry Smith    Input Parameters:
19171713a123SBarry Smith +  ts - the TS context
19181713a123SBarry Smith .  step - current time-step
1919142b95e3SSatish Balay .  ptime - current time
19201713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
19211713a123SBarry Smith 
19221713a123SBarry Smith    Level: intermediate
19231713a123SBarry Smith 
19241713a123SBarry Smith .keywords: TS,  vector, monitor, view
19251713a123SBarry Smith 
19261713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
19271713a123SBarry Smith @*/
1928a7cc72afSBarry Smith PetscErrorCode TSVecViewMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
19291713a123SBarry Smith {
1930dfbe8321SBarry Smith   PetscErrorCode ierr;
19311713a123SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
19321713a123SBarry Smith 
19331713a123SBarry Smith   PetscFunctionBegin;
19341713a123SBarry Smith   if (!viewer) {
19351713a123SBarry Smith     MPI_Comm comm;
19361713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
19371713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
19381713a123SBarry Smith   }
19391713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
19401713a123SBarry Smith   PetscFunctionReturn(0);
19411713a123SBarry Smith }
19421713a123SBarry Smith 
19431713a123SBarry Smith 
19441713a123SBarry Smith 
1945