xref: /petsc/src/ts/interface/ts.c (revision 32077d6df18512e94702c86e7d562481ed0cebd0)
1e090d566SSatish Balay #include "src/ts/tsimpl.h"        /*I "petscts.h"  I*/
2d763cef2SBarry Smith 
3d5ba7fb7SMatthew Knepley /* Logging support */
443c77886SBarry Smith int TS_COOKIE = 0;
543c77886SBarry Smith int 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 */
22bdad233fSMatthew Knepley static int TSSetTypeFromOptions(TS ts)
23bdad233fSMatthew Knepley {
24bdad233fSMatthew Knepley   PetscTruth opt;
252fc52814SBarry Smith   const char *defaultType;
26bdad233fSMatthew Knepley   char       typeName[256];
27bdad233fSMatthew Knepley   int        ierr;
28bdad233fSMatthew Knepley 
29bdad233fSMatthew Knepley   PetscFunctionBegin;
30bdad233fSMatthew Knepley   if (ts->type_name != PETSC_NULL) {
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);
40bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
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 @*/
72bdad233fSMatthew Knepley int TSSetFromOptions(TS ts)
73bdad233fSMatthew Knepley {
74bdad233fSMatthew Knepley   PetscReal  dt;
75bdad233fSMatthew Knepley   PetscTruth opt;
76bdad233fSMatthew Knepley   int        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);
87bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
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);
93bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
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);
97bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
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);
101bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
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 */
109bdad233fSMatthew Knepley   if (ts->ops->setfromoptions != PETSC_NULL) {
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:
11894b7f48cSBarry Smith     if (ts->ksp != PETSC_NULL) {
11994b7f48cSBarry Smith       ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr);
120156fc9a6SMatthew Knepley     }
121bdad233fSMatthew Knepley     break;
122bdad233fSMatthew Knepley   case TS_NONLINEAR:
123156fc9a6SMatthew Knepley     if (ts->snes != PETSC_NULL) {
124bdad233fSMatthew Knepley       ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
125156fc9a6SMatthew Knepley     }
126bdad233fSMatthew Knepley     break;
127bdad233fSMatthew Knepley   default:
128bdad233fSMatthew Knepley     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", ts->problem_type);
129bdad233fSMatthew Knepley   }
130bdad233fSMatthew Knepley 
131bdad233fSMatthew Knepley   PetscFunctionReturn(0);
132bdad233fSMatthew Knepley }
133bdad233fSMatthew Knepley 
134bdad233fSMatthew Knepley #undef  __FUNCT__
135bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
136bdad233fSMatthew Knepley /*@
137bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
138bdad233fSMatthew Knepley 
139bdad233fSMatthew Knepley   Collective on TS
140bdad233fSMatthew Knepley 
141bdad233fSMatthew Knepley   Input Parameter:
142bdad233fSMatthew Knepley . ts - The ts
143bdad233fSMatthew Knepley 
144bdad233fSMatthew Knepley   Level: intermediate
145bdad233fSMatthew Knepley 
146bdad233fSMatthew Knepley .keywords: TS, view, options, database
147bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
148bdad233fSMatthew Knepley @*/
1491836bdbcSSatish Balay int TSViewFromOptions(TS ts,const char title[])
150bdad233fSMatthew Knepley {
151bdad233fSMatthew Knepley   PetscViewer viewer;
152bdad233fSMatthew Knepley   PetscDraw   draw;
153bdad233fSMatthew Knepley   PetscTruth  opt;
154bdad233fSMatthew Knepley   char        typeName[1024];
155e10c49a3SBarry Smith   char        fileName[PETSC_MAX_PATH_LEN];
156bdad233fSMatthew Knepley   int         len;
157bdad233fSMatthew Knepley   int         ierr;
158bdad233fSMatthew Knepley 
159bdad233fSMatthew Knepley   PetscFunctionBegin;
160bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);CHKERRQ(ierr);
161bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
162bdad233fSMatthew Knepley     ierr = PetscOptionsGetString(ts->prefix, "-ts_view", typeName, 1024, &opt);CHKERRQ(ierr);
163bdad233fSMatthew Knepley     ierr = PetscStrlen(typeName, &len);CHKERRQ(ierr);
164bdad233fSMatthew Knepley     if (len > 0) {
165bdad233fSMatthew Knepley       ierr = PetscViewerCreate(ts->comm, &viewer);CHKERRQ(ierr);
166bdad233fSMatthew Knepley       ierr = PetscViewerSetType(viewer, typeName);CHKERRQ(ierr);
167bdad233fSMatthew Knepley       ierr = PetscOptionsGetString(ts->prefix, "-ts_view_file", fileName, 1024, &opt);CHKERRQ(ierr);
168bdad233fSMatthew Knepley       if (opt == PETSC_TRUE) {
169bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, fileName);CHKERRQ(ierr);
170bdad233fSMatthew Knepley       } else {
171bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, ts->name);CHKERRQ(ierr);
172bdad233fSMatthew Knepley       }
173bdad233fSMatthew Knepley       ierr = TSView(ts, viewer);CHKERRQ(ierr);
174bdad233fSMatthew Knepley       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
175bdad233fSMatthew Knepley       ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
176bdad233fSMatthew Knepley     } else {
177bdad233fSMatthew Knepley       ierr = TSView(ts, PETSC_NULL);CHKERRQ(ierr);
178bdad233fSMatthew Knepley     }
179bdad233fSMatthew Knepley   }
180bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);CHKERRQ(ierr);
181bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
182bdad233fSMatthew Knepley     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
183bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
184bdad233fSMatthew Knepley     if (title != PETSC_NULL) {
1851836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
186bdad233fSMatthew Knepley     } else {
187bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject) ts);                                                           CHKERRQ(ierr) ;
1881836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, ts->name);CHKERRQ(ierr);
189bdad233fSMatthew Knepley     }
190bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
191bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
192bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
193bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
194bdad233fSMatthew Knepley   }
195bdad233fSMatthew Knepley   PetscFunctionReturn(0);
196bdad233fSMatthew Knepley }
197bdad233fSMatthew Knepley 
198bdad233fSMatthew Knepley #undef __FUNCT__
1994a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
200a7a1495cSBarry Smith /*@
2018c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
202a7a1495cSBarry Smith       set with TSSetRHSJacobian().
203a7a1495cSBarry Smith 
204a7a1495cSBarry Smith    Collective on TS and Vec
205a7a1495cSBarry Smith 
206a7a1495cSBarry Smith    Input Parameters:
207a7a1495cSBarry Smith +  ts - the SNES context
208a7a1495cSBarry Smith .  t - current timestep
209a7a1495cSBarry Smith -  x - input vector
210a7a1495cSBarry Smith 
211a7a1495cSBarry Smith    Output Parameters:
212a7a1495cSBarry Smith +  A - Jacobian matrix
213a7a1495cSBarry Smith .  B - optional preconditioning matrix
214a7a1495cSBarry Smith -  flag - flag indicating matrix structure
215a7a1495cSBarry Smith 
216a7a1495cSBarry Smith    Notes:
217a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
218a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
219a7a1495cSBarry Smith 
22094b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
221a7a1495cSBarry Smith    flag parameter.
222a7a1495cSBarry Smith 
223a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
224a7a1495cSBarry Smith 
225a7a1495cSBarry Smith    Level: developer
226a7a1495cSBarry Smith 
227a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
228a7a1495cSBarry Smith 
22994b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
230a7a1495cSBarry Smith @*/
23187828ca2SBarry Smith int TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
232a7a1495cSBarry Smith {
233a7a1495cSBarry Smith   int ierr;
234a7a1495cSBarry Smith 
235a7a1495cSBarry Smith   PetscFunctionBegin;
2364482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2374482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
238c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
239a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
24029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
241a7a1495cSBarry Smith   }
242000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
243d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
244a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
245a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
246000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
247a7a1495cSBarry Smith     PetscStackPop;
248d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
249a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2504482741eSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE,4);
2514482741eSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE,5);
252ef66eb69SBarry Smith   } else {
253ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
254ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255ef66eb69SBarry Smith     if (*A != *B) {
256ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258ef66eb69SBarry Smith     }
259ef66eb69SBarry Smith   }
260a7a1495cSBarry Smith   PetscFunctionReturn(0);
261a7a1495cSBarry Smith }
262a7a1495cSBarry Smith 
2634a2ae208SSatish Balay #undef __FUNCT__
2644a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
265d763cef2SBarry Smith /*
266d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
267d763cef2SBarry Smith 
268d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
269d763cef2SBarry Smith    this routine applies the matrix.
270d763cef2SBarry Smith */
27187828ca2SBarry Smith int TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
272d763cef2SBarry Smith {
273d763cef2SBarry Smith   int ierr;
274d763cef2SBarry Smith 
275d763cef2SBarry Smith   PetscFunctionBegin;
2764482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2774482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
2784482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
279d763cef2SBarry Smith 
280d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
281000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
282d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
283000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
284d763cef2SBarry Smith     PetscStackPop;
2857c922b88SBarry Smith   } else {
286000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
287d763cef2SBarry Smith       MatStructure flg;
288d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
289000e7ae3SMatthew Knepley       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
290d763cef2SBarry Smith       PetscStackPop;
291d763cef2SBarry Smith     }
292d763cef2SBarry Smith     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
2937c922b88SBarry Smith   }
294d763cef2SBarry Smith 
295d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
296d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
297d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
298d763cef2SBarry Smith 
299d763cef2SBarry Smith   PetscFunctionReturn(0);
300d763cef2SBarry Smith }
301d763cef2SBarry Smith 
3024a2ae208SSatish Balay #undef __FUNCT__
3034a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
304d763cef2SBarry Smith /*@C
305d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
306d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
307d763cef2SBarry Smith 
308d763cef2SBarry Smith     Collective on TS
309d763cef2SBarry Smith 
310d763cef2SBarry Smith     Input Parameters:
311d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
312d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
313d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
314d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
315d763cef2SBarry Smith 
316d763cef2SBarry Smith     Calling sequence of func:
31787828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
318d763cef2SBarry Smith 
319d763cef2SBarry Smith +   t - current timestep
320d763cef2SBarry Smith .   u - input vector
321d763cef2SBarry Smith .   F - function vector
322d763cef2SBarry Smith -   ctx - [optional] user-defined function context
323d763cef2SBarry Smith 
324d763cef2SBarry Smith     Important:
325d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
326d763cef2SBarry Smith 
327d763cef2SBarry Smith     Level: beginner
328d763cef2SBarry Smith 
329d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
330d763cef2SBarry Smith 
331d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
332d763cef2SBarry Smith @*/
33387828ca2SBarry Smith int TSSetRHSFunction(TS ts,int (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
334d763cef2SBarry Smith {
335d763cef2SBarry Smith   PetscFunctionBegin;
336d763cef2SBarry Smith 
3374482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
338d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
33929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
340d763cef2SBarry Smith   }
341000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
342d763cef2SBarry Smith   ts->funP             = ctx;
343d763cef2SBarry Smith   PetscFunctionReturn(0);
344d763cef2SBarry Smith }
345d763cef2SBarry Smith 
3464a2ae208SSatish Balay #undef __FUNCT__
3474a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
348d763cef2SBarry Smith /*@C
349d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
350d763cef2SBarry Smith    Also sets the location to store A.
351d763cef2SBarry Smith 
352d763cef2SBarry Smith    Collective on TS
353d763cef2SBarry Smith 
354d763cef2SBarry Smith    Input Parameters:
355d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
356d763cef2SBarry Smith .  A   - matrix
357d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
358d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
359d763cef2SBarry Smith          if A is not a function of t.
360d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
361d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
362d763cef2SBarry Smith 
363d763cef2SBarry Smith    Calling sequence of func:
36487828ca2SBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,int *flag,void *ctx);
365d763cef2SBarry Smith 
366d763cef2SBarry Smith +  t - current timestep
367d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
368d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
369d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
37094b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
371d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
372d763cef2SBarry Smith 
373d763cef2SBarry Smith    Notes:
37494b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
375d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
376d763cef2SBarry Smith 
377d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
378d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
379d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
380d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
381d763cef2SBarry Smith    throughout the global iterations.
382d763cef2SBarry Smith 
383d763cef2SBarry Smith    Important:
384d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
385d763cef2SBarry Smith 
386d763cef2SBarry Smith    Level: beginner
387d763cef2SBarry Smith 
388d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
389d763cef2SBarry Smith 
390d763cef2SBarry Smith .seealso: TSSetRHSFunction()
391d763cef2SBarry Smith @*/
39287828ca2SBarry Smith int TSSetRHSMatrix(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
393d763cef2SBarry Smith {
394d763cef2SBarry Smith   PetscFunctionBegin;
3954482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
3964482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
3974482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
398c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
399c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
400d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
40129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
402d763cef2SBarry Smith   }
403d763cef2SBarry Smith 
404000e7ae3SMatthew Knepley   ts->ops->rhsmatrix = f;
405d763cef2SBarry Smith   ts->jacP           = ctx;
406d763cef2SBarry Smith   ts->A              = A;
407d763cef2SBarry Smith   ts->B              = B;
408d763cef2SBarry Smith 
409d763cef2SBarry Smith   PetscFunctionReturn(0);
410d763cef2SBarry Smith }
411d763cef2SBarry Smith 
4124a2ae208SSatish Balay #undef __FUNCT__
4134a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
414d763cef2SBarry Smith /*@C
415d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
416d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
4173c94ec11SBarry Smith    Use TSSetRHSMatrix() for linear problems.
418d763cef2SBarry Smith 
419d763cef2SBarry Smith    Collective on TS
420d763cef2SBarry Smith 
421d763cef2SBarry Smith    Input Parameters:
422d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
423d763cef2SBarry Smith .  A   - Jacobian matrix
424d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
425d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
426d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
427d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
428d763cef2SBarry Smith 
429d763cef2SBarry Smith    Calling sequence of func:
43087828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
431d763cef2SBarry Smith 
432d763cef2SBarry Smith +  t - current timestep
433d763cef2SBarry Smith .  u - input vector
434d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
435d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
436d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
43794b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
438d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
439d763cef2SBarry Smith 
440d763cef2SBarry Smith    Notes:
44194b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
442d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
443d763cef2SBarry Smith 
444d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
445d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
446d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
447d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
448d763cef2SBarry Smith    throughout the global iterations.
449d763cef2SBarry Smith 
450d763cef2SBarry Smith    Level: beginner
451d763cef2SBarry Smith 
452d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
453d763cef2SBarry Smith 
454d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
4553c94ec11SBarry Smith           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetRHSMatrix()
456d763cef2SBarry Smith 
457d763cef2SBarry Smith @*/
45887828ca2SBarry Smith int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
459d763cef2SBarry Smith {
460d763cef2SBarry Smith   PetscFunctionBegin;
4614482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4624482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
4634482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
464c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
465c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
466d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
46729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
468d763cef2SBarry Smith   }
469d763cef2SBarry Smith 
470000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
471d763cef2SBarry Smith   ts->jacP             = ctx;
472d763cef2SBarry Smith   ts->A                = A;
473d763cef2SBarry Smith   ts->B                = B;
474d763cef2SBarry Smith   PetscFunctionReturn(0);
475d763cef2SBarry Smith }
476d763cef2SBarry Smith 
4774a2ae208SSatish Balay #undef __FUNCT__
4784a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
479d763cef2SBarry Smith /*
480d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
481d763cef2SBarry Smith 
482d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
483d763cef2SBarry Smith    this routine applies the matrix.
484d763cef2SBarry Smith */
48587828ca2SBarry Smith int TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
486d763cef2SBarry Smith {
487d763cef2SBarry Smith   int ierr;
488d763cef2SBarry Smith 
489d763cef2SBarry Smith   PetscFunctionBegin;
4904482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4914482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
492c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,x,3);
493d763cef2SBarry Smith 
494000e7ae3SMatthew Knepley   if (ts->ops->rhsbc) {
495d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
496000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
497d763cef2SBarry Smith     PetscStackPop;
498d763cef2SBarry Smith     PetscFunctionReturn(0);
499d763cef2SBarry Smith   }
500d763cef2SBarry Smith 
501d763cef2SBarry Smith   PetscFunctionReturn(0);
502d763cef2SBarry Smith }
503d763cef2SBarry Smith 
5044a2ae208SSatish Balay #undef __FUNCT__
5054a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
506d763cef2SBarry Smith /*@C
507d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
508d763cef2SBarry Smith     boundary conditions for the function F.
509d763cef2SBarry Smith 
510d763cef2SBarry Smith     Collective on TS
511d763cef2SBarry Smith 
512d763cef2SBarry Smith     Input Parameters:
513d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
514d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
515d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
516d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
517d763cef2SBarry Smith 
518d763cef2SBarry Smith     Calling sequence of func:
51987828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
520d763cef2SBarry Smith 
521d763cef2SBarry Smith +   t - current timestep
522d763cef2SBarry Smith .   F - function vector
523d763cef2SBarry Smith -   ctx - [optional] user-defined function context
524d763cef2SBarry Smith 
525d763cef2SBarry Smith     Level: intermediate
526d763cef2SBarry Smith 
527d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
528d763cef2SBarry Smith @*/
52987828ca2SBarry Smith int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,PetscReal,Vec,void*),void *ctx)
530d763cef2SBarry Smith {
531d763cef2SBarry Smith   PetscFunctionBegin;
532d763cef2SBarry Smith 
5334482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
534d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
53529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
536d763cef2SBarry Smith   }
537000e7ae3SMatthew Knepley   ts->ops->rhsbc = f;
538d763cef2SBarry Smith   ts->bcP        = ctx;
539d763cef2SBarry Smith   PetscFunctionReturn(0);
540d763cef2SBarry Smith }
541d763cef2SBarry Smith 
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "TSView"
5447e2c5f70SBarry Smith /*@C
545d763cef2SBarry Smith     TSView - Prints the TS data structure.
546d763cef2SBarry Smith 
5474c49b128SBarry Smith     Collective on TS
548d763cef2SBarry Smith 
549d763cef2SBarry Smith     Input Parameters:
550d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
551d763cef2SBarry Smith -   viewer - visualization context
552d763cef2SBarry Smith 
553d763cef2SBarry Smith     Options Database Key:
554d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
555d763cef2SBarry Smith 
556d763cef2SBarry Smith     Notes:
557d763cef2SBarry Smith     The available visualization contexts include
558b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
559b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
560d763cef2SBarry Smith          output where only the first processor opens
561d763cef2SBarry Smith          the file.  All other processors send their
562d763cef2SBarry Smith          data to the first processor to print.
563d763cef2SBarry Smith 
564d763cef2SBarry Smith     The user can open an alternative visualization context with
565b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
566d763cef2SBarry Smith 
567d763cef2SBarry Smith     Level: beginner
568d763cef2SBarry Smith 
569d763cef2SBarry Smith .keywords: TS, timestep, view
570d763cef2SBarry Smith 
571b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
572d763cef2SBarry Smith @*/
573b0a32e0cSBarry Smith int TSView(TS ts,PetscViewer viewer)
574d763cef2SBarry Smith {
575d763cef2SBarry Smith   int        ierr;
576454a90a3SBarry Smith   char       *type;
577*32077d6dSBarry Smith   PetscTruth iascii,isstring;
578d763cef2SBarry Smith 
579d763cef2SBarry Smith   PetscFunctionBegin;
5804482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
581b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
5824482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
583c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
584fd16b177SBarry Smith 
585*32077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
586b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
587*32077d6dSBarry Smith   if (iascii) {
588b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
589454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
590454a90a3SBarry Smith     if (type) {
591b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
592184914b5SBarry Smith     } else {
593b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
594184914b5SBarry Smith     }
595000e7ae3SMatthew Knepley     if (ts->ops->view) {
596b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
597000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
598b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
599d763cef2SBarry Smith     }
600b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
601b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
602d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
603b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
604d763cef2SBarry Smith     }
605b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
6060f5bd95cSBarry Smith   } else if (isstring) {
607454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
608b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
609d763cef2SBarry Smith   }
610b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
61194b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
612d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
613b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
614d763cef2SBarry Smith   PetscFunctionReturn(0);
615d763cef2SBarry Smith }
616d763cef2SBarry Smith 
617d763cef2SBarry Smith 
6184a2ae208SSatish Balay #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
620d763cef2SBarry Smith /*@C
621d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
622d763cef2SBarry Smith    the timesteppers.
623d763cef2SBarry Smith 
624d763cef2SBarry Smith    Collective on TS
625d763cef2SBarry Smith 
626d763cef2SBarry Smith    Input Parameters:
627d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
628d763cef2SBarry Smith -  usrP - optional user context
629d763cef2SBarry Smith 
630d763cef2SBarry Smith    Level: intermediate
631d763cef2SBarry Smith 
632d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
633d763cef2SBarry Smith 
634d763cef2SBarry Smith .seealso: TSGetApplicationContext()
635d763cef2SBarry Smith @*/
636d763cef2SBarry Smith int TSSetApplicationContext(TS ts,void *usrP)
637d763cef2SBarry Smith {
638d763cef2SBarry Smith   PetscFunctionBegin;
6394482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
640d763cef2SBarry Smith   ts->user = usrP;
641d763cef2SBarry Smith   PetscFunctionReturn(0);
642d763cef2SBarry Smith }
643d763cef2SBarry Smith 
6444a2ae208SSatish Balay #undef __FUNCT__
6454a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
646d763cef2SBarry Smith /*@C
647d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
648d763cef2SBarry Smith     timestepper.
649d763cef2SBarry Smith 
650d763cef2SBarry Smith     Not Collective
651d763cef2SBarry Smith 
652d763cef2SBarry Smith     Input Parameter:
653d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
654d763cef2SBarry Smith 
655d763cef2SBarry Smith     Output Parameter:
656d763cef2SBarry Smith .   usrP - user context
657d763cef2SBarry Smith 
658d763cef2SBarry Smith     Level: intermediate
659d763cef2SBarry Smith 
660d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
661d763cef2SBarry Smith 
662d763cef2SBarry Smith .seealso: TSSetApplicationContext()
663d763cef2SBarry Smith @*/
664d763cef2SBarry Smith int TSGetApplicationContext(TS ts,void **usrP)
665d763cef2SBarry Smith {
666d763cef2SBarry Smith   PetscFunctionBegin;
6674482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
668d763cef2SBarry Smith   *usrP = ts->user;
669d763cef2SBarry Smith   PetscFunctionReturn(0);
670d763cef2SBarry Smith }
671d763cef2SBarry Smith 
6724a2ae208SSatish Balay #undef __FUNCT__
6734a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
674d763cef2SBarry Smith /*@
675d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
676d763cef2SBarry Smith 
677d763cef2SBarry Smith    Not Collective
678d763cef2SBarry Smith 
679d763cef2SBarry Smith    Input Parameter:
680d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
681d763cef2SBarry Smith 
682d763cef2SBarry Smith    Output Parameter:
683d763cef2SBarry Smith .  iter - number steps so far
684d763cef2SBarry Smith 
685d763cef2SBarry Smith    Level: intermediate
686d763cef2SBarry Smith 
687d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
688d763cef2SBarry Smith @*/
689d763cef2SBarry Smith int TSGetTimeStepNumber(TS ts,int* iter)
690d763cef2SBarry Smith {
691d763cef2SBarry Smith   PetscFunctionBegin;
6924482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
6934482741eSBarry Smith   PetscValidIntPointer(iter,2);
694d763cef2SBarry Smith   *iter = ts->steps;
695d763cef2SBarry Smith   PetscFunctionReturn(0);
696d763cef2SBarry Smith }
697d763cef2SBarry Smith 
6984a2ae208SSatish Balay #undef __FUNCT__
6994a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
700d763cef2SBarry Smith /*@
701d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
702d763cef2SBarry Smith    as well as the initial time.
703d763cef2SBarry Smith 
704d763cef2SBarry Smith    Collective on TS
705d763cef2SBarry Smith 
706d763cef2SBarry Smith    Input Parameters:
707d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
708d763cef2SBarry Smith .  initial_time - the initial time
709d763cef2SBarry Smith -  time_step - the size of the timestep
710d763cef2SBarry Smith 
711d763cef2SBarry Smith    Level: intermediate
712d763cef2SBarry Smith 
713d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
714d763cef2SBarry Smith 
715d763cef2SBarry Smith .keywords: TS, set, initial, timestep
716d763cef2SBarry Smith @*/
71787828ca2SBarry Smith int TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
718d763cef2SBarry Smith {
719d763cef2SBarry Smith   PetscFunctionBegin;
7204482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
721d763cef2SBarry Smith   ts->time_step         = time_step;
722d763cef2SBarry Smith   ts->initial_time_step = time_step;
723d763cef2SBarry Smith   ts->ptime             = initial_time;
724d763cef2SBarry Smith   PetscFunctionReturn(0);
725d763cef2SBarry Smith }
726d763cef2SBarry Smith 
7274a2ae208SSatish Balay #undef __FUNCT__
7284a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
729d763cef2SBarry Smith /*@
730d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
731d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
732d763cef2SBarry Smith 
733d763cef2SBarry Smith    Collective on TS
734d763cef2SBarry Smith 
735d763cef2SBarry Smith    Input Parameters:
736d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
737d763cef2SBarry Smith -  time_step - the size of the timestep
738d763cef2SBarry Smith 
739d763cef2SBarry Smith    Level: intermediate
740d763cef2SBarry Smith 
741d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
742d763cef2SBarry Smith 
743d763cef2SBarry Smith .keywords: TS, set, timestep
744d763cef2SBarry Smith @*/
74587828ca2SBarry Smith int TSSetTimeStep(TS ts,PetscReal time_step)
746d763cef2SBarry Smith {
747d763cef2SBarry Smith   PetscFunctionBegin;
7484482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
749d763cef2SBarry Smith   ts->time_step = time_step;
750d763cef2SBarry Smith   PetscFunctionReturn(0);
751d763cef2SBarry Smith }
752d763cef2SBarry Smith 
7534a2ae208SSatish Balay #undef __FUNCT__
7544a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
755d763cef2SBarry Smith /*@
756d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
757d763cef2SBarry Smith 
758d763cef2SBarry Smith    Not Collective
759d763cef2SBarry Smith 
760d763cef2SBarry Smith    Input Parameter:
761d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
762d763cef2SBarry Smith 
763d763cef2SBarry Smith    Output Parameter:
764d763cef2SBarry Smith .  dt - the current timestep size
765d763cef2SBarry Smith 
766d763cef2SBarry Smith    Level: intermediate
767d763cef2SBarry Smith 
768d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
769d763cef2SBarry Smith 
770d763cef2SBarry Smith .keywords: TS, get, timestep
771d763cef2SBarry Smith @*/
77287828ca2SBarry Smith int TSGetTimeStep(TS ts,PetscReal* dt)
773d763cef2SBarry Smith {
774d763cef2SBarry Smith   PetscFunctionBegin;
7754482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
7764482741eSBarry Smith   PetscValidDoublePointer(dt,2);
777d763cef2SBarry Smith   *dt = ts->time_step;
778d763cef2SBarry Smith   PetscFunctionReturn(0);
779d763cef2SBarry Smith }
780d763cef2SBarry Smith 
7814a2ae208SSatish Balay #undef __FUNCT__
7824a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
783d763cef2SBarry Smith /*@C
784d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
785d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
786d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
787d763cef2SBarry Smith    the solution at the next timestep has been calculated.
788d763cef2SBarry Smith 
789d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
790d763cef2SBarry Smith 
791d763cef2SBarry Smith    Input Parameter:
792d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
793d763cef2SBarry Smith 
794d763cef2SBarry Smith    Output Parameter:
795d763cef2SBarry Smith .  v - the vector containing the solution
796d763cef2SBarry Smith 
797d763cef2SBarry Smith    Level: intermediate
798d763cef2SBarry Smith 
799d763cef2SBarry Smith .seealso: TSGetTimeStep()
800d763cef2SBarry Smith 
801d763cef2SBarry Smith .keywords: TS, timestep, get, solution
802d763cef2SBarry Smith @*/
803d763cef2SBarry Smith int TSGetSolution(TS ts,Vec *v)
804d763cef2SBarry Smith {
805d763cef2SBarry Smith   PetscFunctionBegin;
8064482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
8074482741eSBarry Smith   PetscValidPointer(v,2);
808d763cef2SBarry Smith   *v = ts->vec_sol_always;
809d763cef2SBarry Smith   PetscFunctionReturn(0);
810d763cef2SBarry Smith }
811d763cef2SBarry Smith 
812bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
8134a2ae208SSatish Balay #undef __FUNCT__
814bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
815d763cef2SBarry Smith /*@C
816bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
817d763cef2SBarry Smith 
818bdad233fSMatthew Knepley   Not collective
819d763cef2SBarry Smith 
820bdad233fSMatthew Knepley   Input Parameters:
821bdad233fSMatthew Knepley + ts   - The TS
822bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
823d763cef2SBarry Smith .vb
824d763cef2SBarry Smith          U_t = A U
825d763cef2SBarry Smith          U_t = A(t) U
826d763cef2SBarry Smith          U_t = F(t,U)
827d763cef2SBarry Smith .ve
828d763cef2SBarry Smith 
829d763cef2SBarry Smith    Level: beginner
830d763cef2SBarry Smith 
831bdad233fSMatthew Knepley .keywords: TS, problem type
832bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
833d763cef2SBarry Smith @*/
834bdad233fSMatthew Knepley int TSSetProblemType(TS ts, TSProblemType type) {
835d763cef2SBarry Smith   PetscFunctionBegin;
8364482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
837bdad233fSMatthew Knepley   ts->problem_type = type;
838d763cef2SBarry Smith   PetscFunctionReturn(0);
839d763cef2SBarry Smith }
840d763cef2SBarry Smith 
841bdad233fSMatthew Knepley #undef __FUNCT__
842bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
843bdad233fSMatthew Knepley /*@C
844bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
845bdad233fSMatthew Knepley 
846bdad233fSMatthew Knepley   Not collective
847bdad233fSMatthew Knepley 
848bdad233fSMatthew Knepley   Input Parameter:
849bdad233fSMatthew Knepley . ts   - The TS
850bdad233fSMatthew Knepley 
851bdad233fSMatthew Knepley   Output Parameter:
852bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
853bdad233fSMatthew Knepley .vb
854bdad233fSMatthew Knepley          U_t = A U
855bdad233fSMatthew Knepley          U_t = A(t) U
856bdad233fSMatthew Knepley          U_t = F(t,U)
857bdad233fSMatthew Knepley .ve
858bdad233fSMatthew Knepley 
859bdad233fSMatthew Knepley    Level: beginner
860bdad233fSMatthew Knepley 
861bdad233fSMatthew Knepley .keywords: TS, problem type
862bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
863bdad233fSMatthew Knepley @*/
864bdad233fSMatthew Knepley int TSGetProblemType(TS ts, TSProblemType *type) {
865bdad233fSMatthew Knepley   PetscFunctionBegin;
8664482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
8674482741eSBarry Smith   PetscValidIntPointer(type,2);
868bdad233fSMatthew Knepley   *type = ts->problem_type;
869bdad233fSMatthew Knepley   PetscFunctionReturn(0);
870bdad233fSMatthew Knepley }
871d763cef2SBarry Smith 
8724a2ae208SSatish Balay #undef __FUNCT__
8734a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
874d763cef2SBarry Smith /*@
875d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
876d763cef2SBarry Smith    of a timestepper.
877d763cef2SBarry Smith 
878d763cef2SBarry Smith    Collective on TS
879d763cef2SBarry Smith 
880d763cef2SBarry Smith    Input Parameter:
881d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
882d763cef2SBarry Smith 
883d763cef2SBarry Smith    Notes:
884d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
885d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
886d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
887d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
888d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
889d763cef2SBarry Smith 
890d763cef2SBarry Smith    Level: advanced
891d763cef2SBarry Smith 
892d763cef2SBarry Smith .keywords: TS, timestep, setup
893d763cef2SBarry Smith 
894d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
895d763cef2SBarry Smith @*/
896d763cef2SBarry Smith int TSSetUp(TS ts)
897d763cef2SBarry Smith {
898d763cef2SBarry Smith   int ierr;
899d763cef2SBarry Smith 
900d763cef2SBarry Smith   PetscFunctionBegin;
9014482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
90229bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
903d763cef2SBarry Smith   if (!ts->type_name) {
904d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
905d763cef2SBarry Smith   }
906000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
907d763cef2SBarry Smith   ts->setupcalled = 1;
908d763cef2SBarry Smith   PetscFunctionReturn(0);
909d763cef2SBarry Smith }
910d763cef2SBarry Smith 
9114a2ae208SSatish Balay #undef __FUNCT__
9124a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
913d763cef2SBarry Smith /*@C
914d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
915d763cef2SBarry Smith    with TSCreate().
916d763cef2SBarry Smith 
917d763cef2SBarry Smith    Collective on TS
918d763cef2SBarry Smith 
919d763cef2SBarry Smith    Input Parameter:
920d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
921d763cef2SBarry Smith 
922d763cef2SBarry Smith    Level: beginner
923d763cef2SBarry Smith 
924d763cef2SBarry Smith .keywords: TS, timestepper, destroy
925d763cef2SBarry Smith 
926d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
927d763cef2SBarry Smith @*/
928d763cef2SBarry Smith int TSDestroy(TS ts)
929d763cef2SBarry Smith {
930329f5518SBarry Smith   int ierr,i;
931d763cef2SBarry Smith 
932d763cef2SBarry Smith   PetscFunctionBegin;
9334482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
934d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
935d763cef2SBarry Smith 
936be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
9370f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
938be0abb6dSBarry Smith 
93994b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
940d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
941000e7ae3SMatthew Knepley   ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);
942329f5518SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
943329f5518SBarry Smith     if (ts->mdestroy[i]) {
944329f5518SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
945329f5518SBarry Smith     }
946329f5518SBarry Smith   }
947b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)ts);
948d763cef2SBarry Smith   PetscHeaderDestroy((PetscObject)ts);
949d763cef2SBarry Smith   PetscFunctionReturn(0);
950d763cef2SBarry Smith }
951d763cef2SBarry Smith 
9524a2ae208SSatish Balay #undef __FUNCT__
9534a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
954d763cef2SBarry Smith /*@C
955d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
956d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
957d763cef2SBarry Smith 
958d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
959d763cef2SBarry Smith 
960d763cef2SBarry Smith    Input Parameter:
961d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
962d763cef2SBarry Smith 
963d763cef2SBarry Smith    Output Parameter:
964d763cef2SBarry Smith .  snes - the nonlinear solver context
965d763cef2SBarry Smith 
966d763cef2SBarry Smith    Notes:
967d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
968d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
96994b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
970d763cef2SBarry Smith 
971d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
972d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
973d763cef2SBarry Smith 
974d763cef2SBarry Smith    Level: beginner
975d763cef2SBarry Smith 
976d763cef2SBarry Smith .keywords: timestep, get, SNES
977d763cef2SBarry Smith @*/
978d763cef2SBarry Smith int TSGetSNES(TS ts,SNES *snes)
979d763cef2SBarry Smith {
980d763cef2SBarry Smith   PetscFunctionBegin;
9814482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
9824482741eSBarry Smith   PetscValidPointer(snes,2);
98394b7f48cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
984d763cef2SBarry Smith   *snes = ts->snes;
985d763cef2SBarry Smith   PetscFunctionReturn(0);
986d763cef2SBarry Smith }
987d763cef2SBarry Smith 
9884a2ae208SSatish Balay #undef __FUNCT__
98994b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
990d763cef2SBarry Smith /*@C
99194b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
992d763cef2SBarry Smith    a TS (timestepper) context.
993d763cef2SBarry Smith 
99494b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
995d763cef2SBarry Smith 
996d763cef2SBarry Smith    Input Parameter:
997d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
998d763cef2SBarry Smith 
999d763cef2SBarry Smith    Output Parameter:
100094b7f48cSBarry Smith .  ksp - the nonlinear solver context
1001d763cef2SBarry Smith 
1002d763cef2SBarry Smith    Notes:
100394b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1004d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1005d763cef2SBarry Smith    KSP and PC contexts as well.
1006d763cef2SBarry Smith 
100794b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
100894b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1009d763cef2SBarry Smith 
1010d763cef2SBarry Smith    Level: beginner
1011d763cef2SBarry Smith 
101294b7f48cSBarry Smith .keywords: timestep, get, KSP
1013d763cef2SBarry Smith @*/
101494b7f48cSBarry Smith int TSGetKSP(TS ts,KSP *ksp)
1015d763cef2SBarry Smith {
1016d763cef2SBarry Smith   PetscFunctionBegin;
10174482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
10184482741eSBarry Smith   PetscValidPointer(ksp,2);
101929bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
102094b7f48cSBarry Smith   *ksp = ts->ksp;
1021d763cef2SBarry Smith   PetscFunctionReturn(0);
1022d763cef2SBarry Smith }
1023d763cef2SBarry Smith 
1024d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1025d763cef2SBarry Smith 
10264a2ae208SSatish Balay #undef __FUNCT__
1027adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1028adb62b0dSMatthew Knepley /*@
1029adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1030adb62b0dSMatthew Knepley    maximum time for iteration.
1031adb62b0dSMatthew Knepley 
1032adb62b0dSMatthew Knepley    Collective on TS
1033adb62b0dSMatthew Knepley 
1034adb62b0dSMatthew Knepley    Input Parameters:
1035adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1036adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1037adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1038adb62b0dSMatthew Knepley 
1039adb62b0dSMatthew Knepley    Level: intermediate
1040adb62b0dSMatthew Knepley 
1041adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1042adb62b0dSMatthew Knepley @*/
1043adb62b0dSMatthew Knepley int TSGetDuration(TS ts, int *maxsteps, PetscReal *maxtime)
1044adb62b0dSMatthew Knepley {
1045adb62b0dSMatthew Knepley   PetscFunctionBegin;
10464482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1047adb62b0dSMatthew Knepley   if (maxsteps != PETSC_NULL) {
10484482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1049adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1050adb62b0dSMatthew Knepley   }
1051adb62b0dSMatthew Knepley   if (maxtime  != PETSC_NULL) {
10524482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1053adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1054adb62b0dSMatthew Knepley   }
1055adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1056adb62b0dSMatthew Knepley }
1057adb62b0dSMatthew Knepley 
1058adb62b0dSMatthew Knepley #undef __FUNCT__
10594a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1060d763cef2SBarry Smith /*@
1061d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1062d763cef2SBarry Smith    maximum time for iteration.
1063d763cef2SBarry Smith 
1064d763cef2SBarry Smith    Collective on TS
1065d763cef2SBarry Smith 
1066d763cef2SBarry Smith    Input Parameters:
1067d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1068d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1069d763cef2SBarry Smith -  maxtime - final time to iterate to
1070d763cef2SBarry Smith 
1071d763cef2SBarry Smith    Options Database Keys:
1072d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1073d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1074d763cef2SBarry Smith 
1075d763cef2SBarry Smith    Notes:
1076d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1077d763cef2SBarry Smith 
1078d763cef2SBarry Smith    Level: intermediate
1079d763cef2SBarry Smith 
1080d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1081d763cef2SBarry Smith @*/
108287828ca2SBarry Smith int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime)
1083d763cef2SBarry Smith {
1084d763cef2SBarry Smith   PetscFunctionBegin;
10854482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1086d763cef2SBarry Smith   ts->max_steps = maxsteps;
1087d763cef2SBarry Smith   ts->max_time  = maxtime;
1088d763cef2SBarry Smith   PetscFunctionReturn(0);
1089d763cef2SBarry Smith }
1090d763cef2SBarry Smith 
10914a2ae208SSatish Balay #undef __FUNCT__
10924a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1093d763cef2SBarry Smith /*@
1094d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1095d763cef2SBarry Smith    for use by the TS routines.
1096d763cef2SBarry Smith 
1097d763cef2SBarry Smith    Collective on TS and Vec
1098d763cef2SBarry Smith 
1099d763cef2SBarry Smith    Input Parameters:
1100d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1101d763cef2SBarry Smith -  x - the solution vector
1102d763cef2SBarry Smith 
1103d763cef2SBarry Smith    Level: beginner
1104d763cef2SBarry Smith 
1105d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1106d763cef2SBarry Smith @*/
1107d763cef2SBarry Smith int TSSetSolution(TS ts,Vec x)
1108d763cef2SBarry Smith {
1109d763cef2SBarry Smith   PetscFunctionBegin;
11104482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
11114482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1112d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1113d763cef2SBarry Smith   PetscFunctionReturn(0);
1114d763cef2SBarry Smith }
1115d763cef2SBarry Smith 
1116e74ef692SMatthew Knepley #undef __FUNCT__
1117e74ef692SMatthew Knepley #define __FUNCT__ "TSSetRhsBC"
1118ac226902SBarry Smith /*@C
1119000e7ae3SMatthew Knepley   TSSetRhsBC - Sets the function which applies boundary conditions
1120000e7ae3SMatthew Knepley   to the Rhs of each system.
1121000e7ae3SMatthew Knepley 
1122000e7ae3SMatthew Knepley   Collective on TS
1123000e7ae3SMatthew Knepley 
1124000e7ae3SMatthew Knepley   Input Parameters:
1125000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1126000e7ae3SMatthew Knepley - func - The function
1127000e7ae3SMatthew Knepley 
1128000e7ae3SMatthew Knepley   Calling sequence of func:
1129000e7ae3SMatthew Knepley . func (TS ts, Vec rhs, void *ctx);
1130000e7ae3SMatthew Knepley 
1131000e7ae3SMatthew Knepley + rhs - The current rhs vector
1132000e7ae3SMatthew Knepley - ctx - The user-context
1133000e7ae3SMatthew Knepley 
1134000e7ae3SMatthew Knepley   Level: intermediate
1135000e7ae3SMatthew Knepley 
1136000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1137000e7ae3SMatthew Knepley @*/
1138000e7ae3SMatthew Knepley int TSSetRhsBC(TS ts, int (*func)(TS, Vec, void *))
1139000e7ae3SMatthew Knepley {
1140000e7ae3SMatthew Knepley   PetscFunctionBegin;
11414482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1142000e7ae3SMatthew Knepley   ts->ops->applyrhsbc = func;
1143000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1144000e7ae3SMatthew Knepley }
1145000e7ae3SMatthew Knepley 
1146e74ef692SMatthew Knepley #undef __FUNCT__
1147e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultRhsBC"
1148000e7ae3SMatthew Knepley /*@
1149000e7ae3SMatthew Knepley   TSDefaultRhsBC - The default boundary condition function which does nothing.
1150000e7ae3SMatthew Knepley 
1151000e7ae3SMatthew Knepley   Collective on TS
1152000e7ae3SMatthew Knepley 
1153000e7ae3SMatthew Knepley   Input Parameters:
1154000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1155000e7ae3SMatthew Knepley . rhs - The Rhs
1156000e7ae3SMatthew Knepley - ctx - The user-context
1157000e7ae3SMatthew Knepley 
1158000e7ae3SMatthew Knepley   Level: developer
1159000e7ae3SMatthew Knepley 
1160000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1161000e7ae3SMatthew Knepley @*/
1162000e7ae3SMatthew Knepley int TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1163000e7ae3SMatthew Knepley {
1164000e7ae3SMatthew Knepley   PetscFunctionBegin;
1165000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1166000e7ae3SMatthew Knepley }
1167000e7ae3SMatthew Knepley 
1168e74ef692SMatthew Knepley #undef __FUNCT__
1169e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSystemMatrixBC"
1170ac226902SBarry Smith /*@C
1171000e7ae3SMatthew Knepley   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1172000e7ae3SMatthew Knepley   to the system matrix and preconditioner of each system.
1173000e7ae3SMatthew Knepley 
1174000e7ae3SMatthew Knepley   Collective on TS
1175000e7ae3SMatthew Knepley 
1176000e7ae3SMatthew Knepley   Input Parameters:
1177000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1178000e7ae3SMatthew Knepley - func - The function
1179000e7ae3SMatthew Knepley 
1180000e7ae3SMatthew Knepley   Calling sequence of func:
1181000e7ae3SMatthew Knepley . func (TS ts, Mat A, Mat B, void *ctx);
1182000e7ae3SMatthew Knepley 
1183000e7ae3SMatthew Knepley + A   - The current system matrix
1184000e7ae3SMatthew Knepley . B   - The current preconditioner
1185000e7ae3SMatthew Knepley - ctx - The user-context
1186000e7ae3SMatthew Knepley 
1187000e7ae3SMatthew Knepley   Level: intermediate
1188000e7ae3SMatthew Knepley 
1189000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1190000e7ae3SMatthew Knepley @*/
1191000e7ae3SMatthew Knepley int TSSetSystemMatrixBC(TS ts, int (*func)(TS, Mat, Mat, void *))
1192000e7ae3SMatthew Knepley {
1193000e7ae3SMatthew Knepley   PetscFunctionBegin;
11944482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1195000e7ae3SMatthew Knepley   ts->ops->applymatrixbc = func;
1196000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1197000e7ae3SMatthew Knepley }
1198000e7ae3SMatthew Knepley 
1199e74ef692SMatthew Knepley #undef __FUNCT__
1200e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSystemMatrixBC"
1201000e7ae3SMatthew Knepley /*@
1202000e7ae3SMatthew Knepley   TSDefaultSystemMatrixBC - The default boundary condition function which
1203000e7ae3SMatthew Knepley   does nothing.
1204000e7ae3SMatthew Knepley 
1205000e7ae3SMatthew Knepley   Collective on TS
1206000e7ae3SMatthew Knepley 
1207000e7ae3SMatthew Knepley   Input Parameters:
1208000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1209000e7ae3SMatthew Knepley . A   - The system matrix
1210000e7ae3SMatthew Knepley . B   - The preconditioner
1211000e7ae3SMatthew Knepley - ctx - The user-context
1212000e7ae3SMatthew Knepley 
1213000e7ae3SMatthew Knepley   Level: developer
1214000e7ae3SMatthew Knepley 
1215000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1216000e7ae3SMatthew Knepley @*/
1217000e7ae3SMatthew Knepley int TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1218000e7ae3SMatthew Knepley {
1219000e7ae3SMatthew Knepley   PetscFunctionBegin;
1220000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1221000e7ae3SMatthew Knepley }
1222000e7ae3SMatthew Knepley 
1223e74ef692SMatthew Knepley #undef __FUNCT__
1224e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSolutionBC"
1225ac226902SBarry Smith /*@C
1226000e7ae3SMatthew Knepley   TSSetSolutionBC - Sets the function which applies boundary conditions
1227000e7ae3SMatthew Knepley   to the solution of each system. This is necessary in nonlinear systems
1228000e7ae3SMatthew Knepley   which time dependent boundary conditions.
1229000e7ae3SMatthew Knepley 
1230000e7ae3SMatthew Knepley   Collective on TS
1231000e7ae3SMatthew Knepley 
1232000e7ae3SMatthew Knepley   Input Parameters:
1233000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1234000e7ae3SMatthew Knepley - func - The function
1235000e7ae3SMatthew Knepley 
1236000e7ae3SMatthew Knepley   Calling sequence of func:
1237000e7ae3SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
1238000e7ae3SMatthew Knepley 
1239000e7ae3SMatthew Knepley + sol - The current solution vector
1240000e7ae3SMatthew Knepley - ctx - The user-context
1241000e7ae3SMatthew Knepley 
1242000e7ae3SMatthew Knepley   Level: intermediate
1243000e7ae3SMatthew Knepley 
1244000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1245000e7ae3SMatthew Knepley @*/
1246000e7ae3SMatthew Knepley int TSSetSolutionBC(TS ts, int (*func)(TS, Vec, void *))
1247000e7ae3SMatthew Knepley {
1248000e7ae3SMatthew Knepley   PetscFunctionBegin;
12494482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1250000e7ae3SMatthew Knepley   ts->ops->applysolbc = func;
1251000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1252000e7ae3SMatthew Knepley }
1253000e7ae3SMatthew Knepley 
1254e74ef692SMatthew Knepley #undef __FUNCT__
1255e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSolutionBC"
1256000e7ae3SMatthew Knepley /*@
1257000e7ae3SMatthew Knepley   TSDefaultSolutionBC - The default boundary condition function which
1258000e7ae3SMatthew Knepley   does nothing.
1259000e7ae3SMatthew Knepley 
1260000e7ae3SMatthew Knepley   Collective on TS
1261000e7ae3SMatthew Knepley 
1262000e7ae3SMatthew Knepley   Input Parameters:
1263000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1264000e7ae3SMatthew Knepley . sol - The solution
1265000e7ae3SMatthew Knepley - ctx - The user-context
1266000e7ae3SMatthew Knepley 
1267000e7ae3SMatthew Knepley   Level: developer
1268000e7ae3SMatthew Knepley 
1269000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1270000e7ae3SMatthew Knepley @*/
1271000e7ae3SMatthew Knepley int TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1272000e7ae3SMatthew Knepley {
1273000e7ae3SMatthew Knepley   PetscFunctionBegin;
1274000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1275000e7ae3SMatthew Knepley }
1276000e7ae3SMatthew Knepley 
1277e74ef692SMatthew Knepley #undef __FUNCT__
1278e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1279ac226902SBarry Smith /*@C
1280000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1281000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1282000e7ae3SMatthew Knepley 
1283000e7ae3SMatthew Knepley   Collective on TS
1284000e7ae3SMatthew Knepley 
1285000e7ae3SMatthew Knepley   Input Parameters:
1286000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1287000e7ae3SMatthew Knepley - func - The function
1288000e7ae3SMatthew Knepley 
1289000e7ae3SMatthew Knepley   Calling sequence of func:
1290000e7ae3SMatthew Knepley . func (TS ts);
1291000e7ae3SMatthew Knepley 
1292000e7ae3SMatthew Knepley   Level: intermediate
1293000e7ae3SMatthew Knepley 
1294000e7ae3SMatthew Knepley .keywords: TS, timestep
1295000e7ae3SMatthew Knepley @*/
1296000e7ae3SMatthew Knepley int TSSetPreStep(TS ts, int (*func)(TS))
1297000e7ae3SMatthew Knepley {
1298000e7ae3SMatthew Knepley   PetscFunctionBegin;
12994482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1300000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1301000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1302000e7ae3SMatthew Knepley }
1303000e7ae3SMatthew Knepley 
1304e74ef692SMatthew Knepley #undef __FUNCT__
1305e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1306000e7ae3SMatthew Knepley /*@
1307000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1308000e7ae3SMatthew Knepley 
1309000e7ae3SMatthew Knepley   Collective on TS
1310000e7ae3SMatthew Knepley 
1311000e7ae3SMatthew Knepley   Input Parameters:
1312000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1313000e7ae3SMatthew Knepley 
1314000e7ae3SMatthew Knepley   Level: developer
1315000e7ae3SMatthew Knepley 
1316000e7ae3SMatthew Knepley .keywords: TS, timestep
1317000e7ae3SMatthew Knepley @*/
1318000e7ae3SMatthew Knepley int TSDefaultPreStep(TS ts)
1319000e7ae3SMatthew Knepley {
1320000e7ae3SMatthew Knepley   PetscFunctionBegin;
1321000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1322000e7ae3SMatthew Knepley }
1323000e7ae3SMatthew Knepley 
1324e74ef692SMatthew Knepley #undef __FUNCT__
1325e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1326ac226902SBarry Smith /*@C
1327000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1328000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1329000e7ae3SMatthew Knepley   the time step.
1330000e7ae3SMatthew Knepley 
1331000e7ae3SMatthew Knepley   Collective on TS
1332000e7ae3SMatthew Knepley 
1333000e7ae3SMatthew Knepley   Input Parameters:
1334000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1335000e7ae3SMatthew Knepley - func - The function
1336000e7ae3SMatthew Knepley 
1337000e7ae3SMatthew Knepley   Calling sequence of func:
1338000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1339000e7ae3SMatthew Knepley 
1340000e7ae3SMatthew Knepley + t   - The current time
1341000e7ae3SMatthew Knepley - dt  - The current time step
1342000e7ae3SMatthew Knepley 
1343000e7ae3SMatthew Knepley   Level: intermediate
1344000e7ae3SMatthew Knepley 
1345000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1346000e7ae3SMatthew Knepley @*/
134789cef881SMatthew Knepley int TSSetUpdate(TS ts, int (*func)(TS, PetscReal, PetscReal *))
1348000e7ae3SMatthew Knepley {
1349000e7ae3SMatthew Knepley   PetscFunctionBegin;
13504482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1351000e7ae3SMatthew Knepley   ts->ops->update = func;
1352000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1353000e7ae3SMatthew Knepley }
1354000e7ae3SMatthew Knepley 
1355e74ef692SMatthew Knepley #undef __FUNCT__
1356e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1357000e7ae3SMatthew Knepley /*@
1358000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1359000e7ae3SMatthew Knepley 
1360000e7ae3SMatthew Knepley   Collective on TS
1361000e7ae3SMatthew Knepley 
1362000e7ae3SMatthew Knepley   Input Parameters:
1363000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1364000e7ae3SMatthew Knepley - t   - The current time
1365000e7ae3SMatthew Knepley 
1366000e7ae3SMatthew Knepley   Output Parameters:
1367000e7ae3SMatthew Knepley . dt  - The current time step
1368000e7ae3SMatthew Knepley 
1369000e7ae3SMatthew Knepley   Level: developer
1370000e7ae3SMatthew Knepley 
1371000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1372000e7ae3SMatthew Knepley @*/
137389cef881SMatthew Knepley int TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1374000e7ae3SMatthew Knepley {
1375000e7ae3SMatthew Knepley   PetscFunctionBegin;
1376000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1377000e7ae3SMatthew Knepley }
1378000e7ae3SMatthew Knepley 
1379e74ef692SMatthew Knepley #undef __FUNCT__
1380e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1381ac226902SBarry Smith /*@C
1382000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1383000e7ae3SMatthew Knepley   called once at the end of time stepping.
1384000e7ae3SMatthew Knepley 
1385000e7ae3SMatthew Knepley   Collective on TS
1386000e7ae3SMatthew Knepley 
1387000e7ae3SMatthew Knepley   Input Parameters:
1388000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1389000e7ae3SMatthew Knepley - func - The function
1390000e7ae3SMatthew Knepley 
1391000e7ae3SMatthew Knepley   Calling sequence of func:
1392000e7ae3SMatthew Knepley . func (TS ts);
1393000e7ae3SMatthew Knepley 
1394000e7ae3SMatthew Knepley   Level: intermediate
1395000e7ae3SMatthew Knepley 
1396000e7ae3SMatthew Knepley .keywords: TS, timestep
1397000e7ae3SMatthew Knepley @*/
1398000e7ae3SMatthew Knepley int TSSetPostStep(TS ts, int (*func)(TS))
1399000e7ae3SMatthew Knepley {
1400000e7ae3SMatthew Knepley   PetscFunctionBegin;
14014482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1402000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1403000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1404000e7ae3SMatthew Knepley }
1405000e7ae3SMatthew Knepley 
1406e74ef692SMatthew Knepley #undef __FUNCT__
1407e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1408000e7ae3SMatthew Knepley /*@
1409000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1410000e7ae3SMatthew Knepley 
1411000e7ae3SMatthew Knepley   Collective on TS
1412000e7ae3SMatthew Knepley 
1413000e7ae3SMatthew Knepley   Input Parameters:
1414000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1415000e7ae3SMatthew Knepley 
1416000e7ae3SMatthew Knepley   Level: developer
1417000e7ae3SMatthew Knepley 
1418000e7ae3SMatthew Knepley .keywords: TS, timestep
1419000e7ae3SMatthew Knepley @*/
1420000e7ae3SMatthew Knepley int TSDefaultPostStep(TS ts)
1421000e7ae3SMatthew Knepley {
1422000e7ae3SMatthew Knepley   PetscFunctionBegin;
1423000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1424000e7ae3SMatthew Knepley }
1425000e7ae3SMatthew Knepley 
1426d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1427d763cef2SBarry Smith 
14284a2ae208SSatish Balay #undef __FUNCT__
14294a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
1430d763cef2SBarry Smith /*@C
1431d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1432d763cef2SBarry Smith    timestep to display the iteration's  progress.
1433d763cef2SBarry Smith 
1434d763cef2SBarry Smith    Collective on TS
1435d763cef2SBarry Smith 
1436d763cef2SBarry Smith    Input Parameters:
1437d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1438d763cef2SBarry Smith .  func - monitoring routine
1439329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1440b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1441b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1442b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1443d763cef2SBarry Smith 
1444d763cef2SBarry Smith    Calling sequence of func:
144587828ca2SBarry Smith $    int func(TS ts,int steps,PetscReal time,Vec x,void *mctx)
1446d763cef2SBarry Smith 
1447d763cef2SBarry Smith +    ts - the TS context
1448d763cef2SBarry Smith .    steps - iteration number
14491f06c33eSBarry Smith .    time - current time
1450d763cef2SBarry Smith .    x - current iterate
1451d763cef2SBarry Smith -    mctx - [optional] monitoring context
1452d763cef2SBarry Smith 
1453d763cef2SBarry Smith    Notes:
1454d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1455d763cef2SBarry Smith    already has been loaded.
1456d763cef2SBarry Smith 
1457d763cef2SBarry Smith    Level: intermediate
1458d763cef2SBarry Smith 
1459d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1460d763cef2SBarry Smith 
1461d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
1462d763cef2SBarry Smith @*/
146387828ca2SBarry Smith int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*))
1464d763cef2SBarry Smith {
1465d763cef2SBarry Smith   PetscFunctionBegin;
14664482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1467d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
146829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1469d763cef2SBarry Smith   }
1470d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1471329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1472d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1473d763cef2SBarry Smith   PetscFunctionReturn(0);
1474d763cef2SBarry Smith }
1475d763cef2SBarry Smith 
14764a2ae208SSatish Balay #undef __FUNCT__
14774a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
1478d763cef2SBarry Smith /*@C
1479d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1480d763cef2SBarry Smith 
1481d763cef2SBarry Smith    Collective on TS
1482d763cef2SBarry Smith 
1483d763cef2SBarry Smith    Input Parameters:
1484d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1485d763cef2SBarry Smith 
1486d763cef2SBarry Smith    Notes:
1487d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1488d763cef2SBarry Smith 
1489d763cef2SBarry Smith    Level: intermediate
1490d763cef2SBarry Smith 
1491d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1492d763cef2SBarry Smith 
1493d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
1494d763cef2SBarry Smith @*/
1495d763cef2SBarry Smith int TSClearMonitor(TS ts)
1496d763cef2SBarry Smith {
1497d763cef2SBarry Smith   PetscFunctionBegin;
14984482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1499d763cef2SBarry Smith   ts->numbermonitors = 0;
1500d763cef2SBarry Smith   PetscFunctionReturn(0);
1501d763cef2SBarry Smith }
1502d763cef2SBarry Smith 
15034a2ae208SSatish Balay #undef __FUNCT__
15044a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
150587828ca2SBarry Smith int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
1506d763cef2SBarry Smith {
1507d132466eSBarry Smith   int ierr;
1508d132466eSBarry Smith 
1509d763cef2SBarry Smith   PetscFunctionBegin;
1510142b95e3SSatish Balay   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1511d763cef2SBarry Smith   PetscFunctionReturn(0);
1512d763cef2SBarry Smith }
1513d763cef2SBarry Smith 
15144a2ae208SSatish Balay #undef __FUNCT__
15154a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1516d763cef2SBarry Smith /*@
1517d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1518d763cef2SBarry Smith 
1519d763cef2SBarry Smith    Collective on TS
1520d763cef2SBarry Smith 
1521d763cef2SBarry Smith    Input Parameter:
1522d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1523d763cef2SBarry Smith 
1524d763cef2SBarry Smith    Output Parameters:
1525d763cef2SBarry Smith +  steps - number of iterations until termination
1526142b95e3SSatish Balay -  ptime - time until termination
1527d763cef2SBarry Smith 
1528d763cef2SBarry Smith    Level: beginner
1529d763cef2SBarry Smith 
1530d763cef2SBarry Smith .keywords: TS, timestep, solve
1531d763cef2SBarry Smith 
1532d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1533d763cef2SBarry Smith @*/
153487828ca2SBarry Smith int TSStep(TS ts,int *steps,PetscReal *ptime)
1535d763cef2SBarry Smith {
1536f1af5d2fSBarry Smith   int        ierr;
1537d763cef2SBarry Smith 
1538d763cef2SBarry Smith   PetscFunctionBegin;
15394482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1540d405a339SMatthew Knepley   if (!ts->setupcalled) {
1541d405a339SMatthew Knepley     ierr = TSSetUp(ts);CHKERRQ(ierr);
1542d405a339SMatthew Knepley   }
1543d405a339SMatthew Knepley 
1544d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1545d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1546000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1547d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1548d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1549d405a339SMatthew Knepley 
15504bb05414SBarry Smith   if (!PetscPreLoadingOn) {
15514bb05414SBarry Smith     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1552d405a339SMatthew Knepley   }
1553d763cef2SBarry Smith   PetscFunctionReturn(0);
1554d763cef2SBarry Smith }
1555d763cef2SBarry Smith 
15564a2ae208SSatish Balay #undef __FUNCT__
15574a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1558d763cef2SBarry Smith /*
1559d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1560d763cef2SBarry Smith */
156187828ca2SBarry Smith int TSMonitor(TS ts,int step,PetscReal ptime,Vec x)
1562d763cef2SBarry Smith {
1563d763cef2SBarry Smith   int i,ierr,n = ts->numbermonitors;
1564d763cef2SBarry Smith 
1565d763cef2SBarry Smith   PetscFunctionBegin;
1566d763cef2SBarry Smith   for (i=0; i<n; i++) {
1567142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1568d763cef2SBarry Smith   }
1569d763cef2SBarry Smith   PetscFunctionReturn(0);
1570d763cef2SBarry Smith }
1571d763cef2SBarry Smith 
1572d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1573d763cef2SBarry Smith 
15744a2ae208SSatish Balay #undef __FUNCT__
15754a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1576d763cef2SBarry Smith /*@C
1577d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1578d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1579d763cef2SBarry Smith 
1580d763cef2SBarry Smith    Collective on TS
1581d763cef2SBarry Smith 
1582d763cef2SBarry Smith    Input Parameters:
1583d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1584d763cef2SBarry Smith .  label - the title to put in the title bar
15857c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1586d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1587d763cef2SBarry Smith 
1588d763cef2SBarry Smith    Output Parameter:
1589d763cef2SBarry Smith .  draw - the drawing context
1590d763cef2SBarry Smith 
1591d763cef2SBarry Smith    Options Database Key:
1592d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1593d763cef2SBarry Smith 
1594d763cef2SBarry Smith    Notes:
1595b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1596d763cef2SBarry Smith 
1597d763cef2SBarry Smith    Level: intermediate
1598d763cef2SBarry Smith 
15997c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1600d763cef2SBarry Smith 
1601d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
16027c922b88SBarry Smith 
1603d763cef2SBarry Smith @*/
16041836bdbcSSatish Balay int TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1605d763cef2SBarry Smith {
1606b0a32e0cSBarry Smith   PetscDraw win;
1607d763cef2SBarry Smith   int       ierr;
1608d763cef2SBarry Smith 
1609d763cef2SBarry Smith   PetscFunctionBegin;
1610b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1611b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1612b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1613b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1614d763cef2SBarry Smith 
1615b0a32e0cSBarry Smith   PetscLogObjectParent(*draw,win);
1616d763cef2SBarry Smith   PetscFunctionReturn(0);
1617d763cef2SBarry Smith }
1618d763cef2SBarry Smith 
16194a2ae208SSatish Balay #undef __FUNCT__
16204a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
162187828ca2SBarry Smith int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx)
1622d763cef2SBarry Smith {
1623b0a32e0cSBarry Smith   PetscDrawLG lg = (PetscDrawLG) monctx;
162487828ca2SBarry Smith   PetscReal      x,y = ptime;
1625d763cef2SBarry Smith   int         ierr;
1626d763cef2SBarry Smith 
1627d763cef2SBarry Smith   PetscFunctionBegin;
16287c922b88SBarry Smith   if (!monctx) {
16297c922b88SBarry Smith     MPI_Comm    comm;
1630b0a32e0cSBarry Smith     PetscViewer viewer;
16317c922b88SBarry Smith 
16327c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1633b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1634b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
16357c922b88SBarry Smith   }
16367c922b88SBarry Smith 
1637b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
163887828ca2SBarry Smith   x = (PetscReal)n;
1639b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1640d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1641b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1642d763cef2SBarry Smith   }
1643d763cef2SBarry Smith   PetscFunctionReturn(0);
1644d763cef2SBarry Smith }
1645d763cef2SBarry Smith 
16464a2ae208SSatish Balay #undef __FUNCT__
16474a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1648d763cef2SBarry Smith /*@C
1649d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1650d763cef2SBarry Smith    with TSLGMonitorCreate().
1651d763cef2SBarry Smith 
1652b0a32e0cSBarry Smith    Collective on PetscDrawLG
1653d763cef2SBarry Smith 
1654d763cef2SBarry Smith    Input Parameter:
1655d763cef2SBarry Smith .  draw - the drawing context
1656d763cef2SBarry Smith 
1657d763cef2SBarry Smith    Level: intermediate
1658d763cef2SBarry Smith 
1659d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1660d763cef2SBarry Smith 
1661d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1662d763cef2SBarry Smith @*/
1663b0a32e0cSBarry Smith int TSLGMonitorDestroy(PetscDrawLG drawlg)
1664d763cef2SBarry Smith {
1665b0a32e0cSBarry Smith   PetscDraw draw;
1666d763cef2SBarry Smith   int       ierr;
1667d763cef2SBarry Smith 
1668d763cef2SBarry Smith   PetscFunctionBegin;
1669b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1670b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1671b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1672d763cef2SBarry Smith   PetscFunctionReturn(0);
1673d763cef2SBarry Smith }
1674d763cef2SBarry Smith 
16754a2ae208SSatish Balay #undef __FUNCT__
16764a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1677d763cef2SBarry Smith /*@
1678d763cef2SBarry Smith    TSGetTime - Gets the current time.
1679d763cef2SBarry Smith 
1680d763cef2SBarry Smith    Not Collective
1681d763cef2SBarry Smith 
1682d763cef2SBarry Smith    Input Parameter:
1683d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1684d763cef2SBarry Smith 
1685d763cef2SBarry Smith    Output Parameter:
1686d763cef2SBarry Smith .  t  - the current time
1687d763cef2SBarry Smith 
1688d763cef2SBarry Smith    Contributed by: Matthew Knepley
1689d763cef2SBarry Smith 
1690d763cef2SBarry Smith    Level: beginner
1691d763cef2SBarry Smith 
1692d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1693d763cef2SBarry Smith 
1694d763cef2SBarry Smith .keywords: TS, get, time
1695d763cef2SBarry Smith @*/
169687828ca2SBarry Smith int TSGetTime(TS ts,PetscReal* t)
1697d763cef2SBarry Smith {
1698d763cef2SBarry Smith   PetscFunctionBegin;
16994482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
17004482741eSBarry Smith   PetscValidDoublePointer(t,2);
1701d763cef2SBarry Smith   *t = ts->ptime;
1702d763cef2SBarry Smith   PetscFunctionReturn(0);
1703d763cef2SBarry Smith }
1704d763cef2SBarry Smith 
17054a2ae208SSatish Balay #undef __FUNCT__
17064a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1707d763cef2SBarry Smith /*@C
1708d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1709d763cef2SBarry Smith    TS options in the database.
1710d763cef2SBarry Smith 
1711d763cef2SBarry Smith    Collective on TS
1712d763cef2SBarry Smith 
1713d763cef2SBarry Smith    Input Parameter:
1714d763cef2SBarry Smith +  ts     - The TS context
1715d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1716d763cef2SBarry Smith 
1717d763cef2SBarry Smith    Notes:
1718d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1719d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1720d763cef2SBarry Smith    hyphen.
1721d763cef2SBarry Smith 
1722d763cef2SBarry Smith    Contributed by: Matthew Knepley
1723d763cef2SBarry Smith 
1724d763cef2SBarry Smith    Level: advanced
1725d763cef2SBarry Smith 
1726d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1727d763cef2SBarry Smith 
1728d763cef2SBarry Smith .seealso: TSSetFromOptions()
1729d763cef2SBarry Smith 
1730d763cef2SBarry Smith @*/
17311836bdbcSSatish Balay int TSSetOptionsPrefix(TS ts,const char prefix[])
1732d763cef2SBarry Smith {
1733d763cef2SBarry Smith   int ierr;
1734d763cef2SBarry Smith 
1735d763cef2SBarry Smith   PetscFunctionBegin;
17364482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1737d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1738d763cef2SBarry Smith   switch(ts->problem_type) {
1739d763cef2SBarry Smith     case TS_NONLINEAR:
1740d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1741d763cef2SBarry Smith       break;
1742d763cef2SBarry Smith     case TS_LINEAR:
174394b7f48cSBarry Smith       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1744d763cef2SBarry Smith       break;
1745d763cef2SBarry Smith   }
1746d763cef2SBarry Smith   PetscFunctionReturn(0);
1747d763cef2SBarry Smith }
1748d763cef2SBarry Smith 
1749d763cef2SBarry Smith 
17504a2ae208SSatish Balay #undef __FUNCT__
17514a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1752d763cef2SBarry Smith /*@C
1753d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1754d763cef2SBarry Smith    TS options in the database.
1755d763cef2SBarry Smith 
1756d763cef2SBarry Smith    Collective on TS
1757d763cef2SBarry Smith 
1758d763cef2SBarry Smith    Input Parameter:
1759d763cef2SBarry Smith +  ts     - The TS context
1760d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1761d763cef2SBarry Smith 
1762d763cef2SBarry Smith    Notes:
1763d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1764d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1765d763cef2SBarry Smith    hyphen.
1766d763cef2SBarry Smith 
1767d763cef2SBarry Smith    Contributed by: Matthew Knepley
1768d763cef2SBarry Smith 
1769d763cef2SBarry Smith    Level: advanced
1770d763cef2SBarry Smith 
1771d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1772d763cef2SBarry Smith 
1773d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1774d763cef2SBarry Smith 
1775d763cef2SBarry Smith @*/
17761836bdbcSSatish Balay int TSAppendOptionsPrefix(TS ts,const char prefix[])
1777d763cef2SBarry Smith {
1778d763cef2SBarry Smith   int ierr;
1779d763cef2SBarry Smith 
1780d763cef2SBarry Smith   PetscFunctionBegin;
17814482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1782d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1783d763cef2SBarry Smith   switch(ts->problem_type) {
1784d763cef2SBarry Smith     case TS_NONLINEAR:
1785d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1786d763cef2SBarry Smith       break;
1787d763cef2SBarry Smith     case TS_LINEAR:
178894b7f48cSBarry Smith       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1789d763cef2SBarry Smith       break;
1790d763cef2SBarry Smith   }
1791d763cef2SBarry Smith   PetscFunctionReturn(0);
1792d763cef2SBarry Smith }
1793d763cef2SBarry Smith 
17944a2ae208SSatish Balay #undef __FUNCT__
17954a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1796d763cef2SBarry Smith /*@C
1797d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1798d763cef2SBarry Smith    TS options in the database.
1799d763cef2SBarry Smith 
1800d763cef2SBarry Smith    Not Collective
1801d763cef2SBarry Smith 
1802d763cef2SBarry Smith    Input Parameter:
1803d763cef2SBarry Smith .  ts - The TS context
1804d763cef2SBarry Smith 
1805d763cef2SBarry Smith    Output Parameter:
1806d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1807d763cef2SBarry Smith 
1808d763cef2SBarry Smith    Contributed by: Matthew Knepley
1809d763cef2SBarry Smith 
1810d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1811d763cef2SBarry Smith    sufficient length to hold the prefix.
1812d763cef2SBarry Smith 
1813d763cef2SBarry Smith    Level: intermediate
1814d763cef2SBarry Smith 
1815d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1816d763cef2SBarry Smith 
1817d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1818d763cef2SBarry Smith @*/
18191836bdbcSSatish Balay int TSGetOptionsPrefix(TS ts,char *prefix[])
1820d763cef2SBarry Smith {
1821d763cef2SBarry Smith   int ierr;
1822d763cef2SBarry Smith 
1823d763cef2SBarry Smith   PetscFunctionBegin;
18244482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
18254482741eSBarry Smith   PetscValidPointer(prefix,2);
1826d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1827d763cef2SBarry Smith   PetscFunctionReturn(0);
1828d763cef2SBarry Smith }
1829d763cef2SBarry Smith 
18304a2ae208SSatish Balay #undef __FUNCT__
18314a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1832d763cef2SBarry Smith /*@C
1833d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1834d763cef2SBarry Smith 
1835d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1836d763cef2SBarry Smith 
1837d763cef2SBarry Smith    Input Parameter:
1838d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1839d763cef2SBarry Smith 
1840d763cef2SBarry Smith    Output Parameters:
1841d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1842d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1843d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1844d763cef2SBarry Smith 
1845d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1846d763cef2SBarry Smith 
1847d763cef2SBarry Smith    Contributed by: Matthew Knepley
1848d763cef2SBarry Smith 
1849d763cef2SBarry Smith    Level: intermediate
1850d763cef2SBarry Smith 
1851d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1852d763cef2SBarry Smith 
1853d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1854d763cef2SBarry Smith 
1855d763cef2SBarry Smith @*/
1856d763cef2SBarry Smith int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1857d763cef2SBarry Smith {
1858d763cef2SBarry Smith   PetscFunctionBegin;
18594482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1860d763cef2SBarry Smith   if (A)   *A = ts->A;
1861d763cef2SBarry Smith   if (M)   *M = ts->B;
1862d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1863d763cef2SBarry Smith   PetscFunctionReturn(0);
1864d763cef2SBarry Smith }
1865d763cef2SBarry Smith 
18664a2ae208SSatish Balay #undef __FUNCT__
18674a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1868d763cef2SBarry Smith /*@C
1869d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1870d763cef2SBarry Smith 
1871d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1872d763cef2SBarry Smith 
1873d763cef2SBarry Smith    Input Parameter:
1874d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1875d763cef2SBarry Smith 
1876d763cef2SBarry Smith    Output Parameters:
1877d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1878d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1879d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1880d763cef2SBarry Smith 
1881d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1882d763cef2SBarry Smith 
1883d763cef2SBarry Smith    Contributed by: Matthew Knepley
1884d763cef2SBarry Smith 
1885d763cef2SBarry Smith    Level: intermediate
1886d763cef2SBarry Smith 
1887d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1888d763cef2SBarry Smith 
1889d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1890d763cef2SBarry Smith @*/
1891d763cef2SBarry Smith int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1892d763cef2SBarry Smith {
1893d763cef2SBarry Smith   int ierr;
1894d763cef2SBarry Smith 
1895d763cef2SBarry Smith   PetscFunctionBegin;
1896d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1897d763cef2SBarry Smith   PetscFunctionReturn(0);
1898d763cef2SBarry Smith }
1899d763cef2SBarry Smith 
19001713a123SBarry Smith #undef __FUNCT__
19011713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
19021713a123SBarry Smith /*@C
19031713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
19041713a123SBarry Smith    VecView() for the solution at each timestep
19051713a123SBarry Smith 
19061713a123SBarry Smith    Collective on TS
19071713a123SBarry Smith 
19081713a123SBarry Smith    Input Parameters:
19091713a123SBarry Smith +  ts - the TS context
19101713a123SBarry Smith .  step - current time-step
1911142b95e3SSatish Balay .  ptime - current time
19121713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
19131713a123SBarry Smith 
19141713a123SBarry Smith    Level: intermediate
19151713a123SBarry Smith 
19161713a123SBarry Smith .keywords: TS,  vector, monitor, view
19171713a123SBarry Smith 
19181713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
19191713a123SBarry Smith @*/
1920142b95e3SSatish Balay int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
19211713a123SBarry Smith {
19221713a123SBarry Smith   int         ierr;
19231713a123SBarry Smith   PetscViewer viewer = (PetscViewer) dummy;
19241713a123SBarry Smith 
19251713a123SBarry Smith   PetscFunctionBegin;
19261713a123SBarry Smith   if (!viewer) {
19271713a123SBarry Smith     MPI_Comm comm;
19281713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
19291713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
19301713a123SBarry Smith   }
19311713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
19321713a123SBarry Smith   PetscFunctionReturn(0);
19331713a123SBarry Smith }
19341713a123SBarry Smith 
19351713a123SBarry Smith 
19361713a123SBarry Smith 
1937