xref: /petsc/src/ts/interface/ts.c (revision 156fc9a69ca21beefd93f4a5b0c292278ce83998)
173f4d377SMatthew Knepley /* $Id: ts.c,v 1.43 2001/09/07 20:12:01 bsmith Exp $ */
2e090d566SSatish Balay #include "src/ts/tsimpl.h"        /*I "petscts.h"  I*/
3d763cef2SBarry Smith 
4d5ba7fb7SMatthew Knepley /* Logging support */
58ba1e511SMatthew Knepley int TS_COOKIE;
6d5ba7fb7SMatthew Knepley int TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
7d405a339SMatthew Knepley 
84a2ae208SSatish Balay #undef __FUNCT__
9bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
10bdad233fSMatthew Knepley /*
11bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
12bdad233fSMatthew Knepley 
13bdad233fSMatthew Knepley   Collective on TS
14bdad233fSMatthew Knepley 
15bdad233fSMatthew Knepley   Input Parameter:
16bdad233fSMatthew Knepley . ts - The ts
17bdad233fSMatthew Knepley 
18bdad233fSMatthew Knepley   Level: intermediate
19bdad233fSMatthew Knepley 
20bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
21bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
22bdad233fSMatthew Knepley */
23bdad233fSMatthew Knepley static int TSSetTypeFromOptions(TS ts)
24bdad233fSMatthew Knepley {
25bdad233fSMatthew Knepley   PetscTruth opt;
26bdad233fSMatthew Knepley   char      *defaultType;
27bdad233fSMatthew Knepley   char       typeName[256];
28bdad233fSMatthew Knepley   int        ierr;
29bdad233fSMatthew Knepley 
30bdad233fSMatthew Knepley   PetscFunctionBegin;
31bdad233fSMatthew Knepley   if (ts->type_name != PETSC_NULL) {
32bdad233fSMatthew Knepley     defaultType = ts->type_name;
33bdad233fSMatthew Knepley   } else {
34bdad233fSMatthew Knepley     defaultType = TS_EULER;
35bdad233fSMatthew Knepley   }
36bdad233fSMatthew Knepley 
37bdad233fSMatthew Knepley   if (!TSRegisterAllCalled) {
38bdad233fSMatthew Knepley     ierr = TSRegisterAll(PETSC_NULL);                                                                     CHKERRQ(ierr);
39bdad233fSMatthew Knepley   }
40bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
41bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
42bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);                                                                       CHKERRQ(ierr);
43bdad233fSMatthew Knepley   } else {
44bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);                                                                    CHKERRQ(ierr);
45bdad233fSMatthew Knepley   }
46bdad233fSMatthew Knepley   PetscFunctionReturn(0);
47bdad233fSMatthew Knepley }
48bdad233fSMatthew Knepley 
49bdad233fSMatthew Knepley #undef __FUNCT__
50bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
51bdad233fSMatthew Knepley /*@
52bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
53bdad233fSMatthew Knepley 
54bdad233fSMatthew Knepley    Collective on TS
55bdad233fSMatthew Knepley 
56bdad233fSMatthew Knepley    Input Parameter:
57bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
58bdad233fSMatthew Knepley 
59bdad233fSMatthew Knepley    Options Database Keys:
60bdad233fSMatthew Knepley +  -ts_type <type> - TS_EULER, TS_BEULER, TS_PVODE, TS_PSEUDO, TS_CRANK_NICHOLSON
61bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
62bdad233fSMatthew Knepley .  -ts_max_time time - maximum time to compute to
63bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
64bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
65bdad233fSMatthew Knepley -  -ts_xmonitor - plot information at each timestep
66bdad233fSMatthew Knepley 
67bdad233fSMatthew Knepley    Level: beginner
68bdad233fSMatthew Knepley 
69bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
70bdad233fSMatthew Knepley 
71bdad233fSMatthew Knepley .seealso: TSGetType
72bdad233fSMatthew Knepley @*/
73bdad233fSMatthew Knepley int TSSetFromOptions(TS ts)
74bdad233fSMatthew Knepley {
75bdad233fSMatthew Knepley   PetscReal  dt;
76bdad233fSMatthew Knepley   PetscTruth opt;
77bdad233fSMatthew Knepley   int        ierr;
78bdad233fSMatthew Knepley 
79bdad233fSMatthew Knepley   PetscFunctionBegin;
80bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
81bdad233fSMatthew Knepley   ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");                              CHKERRQ(ierr);
82bdad233fSMatthew Knepley 
83bdad233fSMatthew Knepley   /* Handle generic TS options */
84bdad233fSMatthew Knepley   ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
85bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
86bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
87bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
88bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
89bdad233fSMatthew Knepley     ts->initial_time_step = ts->time_step = dt;
90bdad233fSMatthew Knepley   }
91bdad233fSMatthew Knepley 
92bdad233fSMatthew Knepley   /* Monitor options */
93bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_monitor","Monitor timestep size","TSDefaultMonitor",&opt);               CHKERRQ(ierr);
94bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
95bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSDefaultMonitor,PETSC_NULL,PETSC_NULL);                                     CHKERRQ(ierr);
96bdad233fSMatthew Knepley     }
97bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_xmonitor","Monitor timestep size graphically","TSLGMonitor",&opt);       CHKERRQ(ierr);
98bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
99bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSLGMonitor,PETSC_NULL,PETSC_NULL);                                          CHKERRQ(ierr);
100bdad233fSMatthew Knepley     }
101bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_vecmonitor","Monitor solution graphically","TSVecViewMonitor",&opt);     CHKERRQ(ierr);
102bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
103bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSVecViewMonitor,PETSC_NULL,PETSC_NULL);                                     CHKERRQ(ierr);
104bdad233fSMatthew Knepley     }
105bdad233fSMatthew Knepley 
106bdad233fSMatthew Knepley   /* Handle TS type options */
107bdad233fSMatthew Knepley   ierr = TSSetTypeFromOptions(ts);                                                                        CHKERRQ(ierr);
108bdad233fSMatthew Knepley 
109bdad233fSMatthew Knepley   /* Handle specific TS options */
110bdad233fSMatthew Knepley   if (ts->ops->setfromoptions != PETSC_NULL) {
111bdad233fSMatthew Knepley     ierr = (*ts->ops->setfromoptions)(ts);                                                                CHKERRQ(ierr);
112bdad233fSMatthew Knepley   }
113bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();                                                                               CHKERRQ(ierr);
114bdad233fSMatthew Knepley 
115bdad233fSMatthew Knepley   /* Handle subobject options */
116bdad233fSMatthew Knepley   switch(ts->problem_type) {
117*156fc9a6SMatthew Knepley     /* Should check for implicit/explicit */
118bdad233fSMatthew Knepley   case TS_LINEAR:
119*156fc9a6SMatthew Knepley     if (ts->sles != PETSC_NULL) {
120bdad233fSMatthew Knepley       ierr = SLESSetFromOptions(ts->sles);                                                                CHKERRQ(ierr);
121*156fc9a6SMatthew Knepley     }
122bdad233fSMatthew Knepley     break;
123bdad233fSMatthew Knepley   case TS_NONLINEAR:
124*156fc9a6SMatthew Knepley     if (ts->snes != PETSC_NULL) {
125bdad233fSMatthew Knepley       ierr = SNESSetFromOptions(ts->snes);                                                                CHKERRQ(ierr);
126*156fc9a6SMatthew Knepley     }
127bdad233fSMatthew Knepley     break;
128bdad233fSMatthew Knepley   default:
129bdad233fSMatthew Knepley     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", ts->problem_type);
130bdad233fSMatthew Knepley   }
131bdad233fSMatthew Knepley 
132bdad233fSMatthew Knepley   ierr = TSViewFromOptions(ts, ts->name);                                                                 CHKERRQ(ierr);
133bdad233fSMatthew Knepley   PetscFunctionReturn(0);
134bdad233fSMatthew Knepley }
135bdad233fSMatthew Knepley 
136bdad233fSMatthew Knepley #undef  __FUNCT__
137bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
138bdad233fSMatthew Knepley /*@
139bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
140bdad233fSMatthew Knepley 
141bdad233fSMatthew Knepley   Collective on TS
142bdad233fSMatthew Knepley 
143bdad233fSMatthew Knepley   Input Parameter:
144bdad233fSMatthew Knepley . ts - The ts
145bdad233fSMatthew Knepley 
146bdad233fSMatthew Knepley   Level: intermediate
147bdad233fSMatthew Knepley 
148bdad233fSMatthew Knepley .keywords: TS, view, options, database
149bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
150bdad233fSMatthew Knepley @*/
151bdad233fSMatthew Knepley int TSViewFromOptions(TS ts, char *title)
152bdad233fSMatthew Knepley {
153bdad233fSMatthew Knepley   PetscViewer viewer;
154bdad233fSMatthew Knepley   PetscDraw   draw;
155bdad233fSMatthew Knepley   PetscTruth  opt;
156bdad233fSMatthew Knepley   char       *titleStr;
157bdad233fSMatthew Knepley   char        typeName[1024];
158bdad233fSMatthew Knepley   char        fileName[1024];
159bdad233fSMatthew Knepley   int         len;
160bdad233fSMatthew Knepley   int         ierr;
161bdad233fSMatthew Knepley 
162bdad233fSMatthew Knepley   PetscFunctionBegin;
163bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);                                               CHKERRQ(ierr);
164bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
165bdad233fSMatthew Knepley     ierr = PetscOptionsGetString(ts->prefix, "-ts_view", typeName, 1024, &opt);                           CHKERRQ(ierr);
166bdad233fSMatthew Knepley     ierr = PetscStrlen(typeName, &len);                                                                   CHKERRQ(ierr);
167bdad233fSMatthew Knepley     if (len > 0) {
168bdad233fSMatthew Knepley       ierr = PetscViewerCreate(ts->comm, &viewer);                                                        CHKERRQ(ierr);
169bdad233fSMatthew Knepley       ierr = PetscViewerSetType(viewer, typeName);                                                        CHKERRQ(ierr);
170bdad233fSMatthew Knepley       ierr = PetscOptionsGetString(ts->prefix, "-ts_view_file", fileName, 1024, &opt);                    CHKERRQ(ierr);
171bdad233fSMatthew Knepley       if (opt == PETSC_TRUE) {
172bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, fileName);                                                  CHKERRQ(ierr);
173bdad233fSMatthew Knepley       } else {
174bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, ts->name);                                                  CHKERRQ(ierr);
175bdad233fSMatthew Knepley       }
176bdad233fSMatthew Knepley       ierr = TSView(ts, viewer);                                                                          CHKERRQ(ierr);
177bdad233fSMatthew Knepley       ierr = PetscViewerFlush(viewer);                                                                    CHKERRQ(ierr);
178bdad233fSMatthew Knepley       ierr = PetscViewerDestroy(viewer);                                                                  CHKERRQ(ierr);
179bdad233fSMatthew Knepley     } else {
180bdad233fSMatthew Knepley       ierr = TSView(ts, PETSC_NULL);                                                                      CHKERRQ(ierr);
181bdad233fSMatthew Knepley     }
182bdad233fSMatthew Knepley   }
183bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);                                          CHKERRQ(ierr);
184bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
185bdad233fSMatthew Knepley     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);                                  CHKERRQ(ierr);
186bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);                                                      CHKERRQ(ierr);
187bdad233fSMatthew Knepley     if (title != PETSC_NULL) {
188bdad233fSMatthew Knepley       titleStr = title;
189bdad233fSMatthew Knepley     } else {
190bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject) ts);                                                           CHKERRQ(ierr) ;
191bdad233fSMatthew Knepley       titleStr = ts->name;
192bdad233fSMatthew Knepley     }
193bdad233fSMatthew Knepley     ierr = PetscDrawSetTitle(draw, titleStr);                                                             CHKERRQ(ierr);
194bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);                                                                            CHKERRQ(ierr);
195bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);                                                                      CHKERRQ(ierr);
196bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);                                                                          CHKERRQ(ierr);
197bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);                                                                    CHKERRQ(ierr);
198bdad233fSMatthew Knepley   }
199bdad233fSMatthew Knepley   PetscFunctionReturn(0);
200bdad233fSMatthew Knepley }
201bdad233fSMatthew Knepley 
202bdad233fSMatthew Knepley #undef __FUNCT__
2034a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
204a7a1495cSBarry Smith /*@
2058c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
206a7a1495cSBarry Smith       set with TSSetRHSJacobian().
207a7a1495cSBarry Smith 
208a7a1495cSBarry Smith    Collective on TS and Vec
209a7a1495cSBarry Smith 
210a7a1495cSBarry Smith    Input Parameters:
211a7a1495cSBarry Smith +  ts - the SNES context
212a7a1495cSBarry Smith .  t - current timestep
213a7a1495cSBarry Smith -  x - input vector
214a7a1495cSBarry Smith 
215a7a1495cSBarry Smith    Output Parameters:
216a7a1495cSBarry Smith +  A - Jacobian matrix
217a7a1495cSBarry Smith .  B - optional preconditioning matrix
218a7a1495cSBarry Smith -  flag - flag indicating matrix structure
219a7a1495cSBarry Smith 
220a7a1495cSBarry Smith    Notes:
221a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
222a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
223a7a1495cSBarry Smith 
224a7a1495cSBarry Smith    See SLESSetOperators() for important information about setting the
225a7a1495cSBarry Smith    flag parameter.
226a7a1495cSBarry Smith 
227a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
228a7a1495cSBarry Smith 
229a7a1495cSBarry Smith    Level: developer
230a7a1495cSBarry Smith 
231a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
232a7a1495cSBarry Smith 
233a7a1495cSBarry Smith .seealso:  TSSetRHSJacobian(), SLESSetOperators()
234a7a1495cSBarry Smith @*/
23587828ca2SBarry Smith int TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
236a7a1495cSBarry Smith {
237a7a1495cSBarry Smith   int ierr;
238a7a1495cSBarry Smith 
239a7a1495cSBarry Smith   PetscFunctionBegin;
240a7a1495cSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
241a7a1495cSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE);
242a7a1495cSBarry Smith   PetscCheckSameComm(ts,X);
243a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
24429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
245a7a1495cSBarry Smith   }
246000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
247d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
248a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
249a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
250000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
251a7a1495cSBarry Smith     PetscStackPop;
252d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
253a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
254a7a1495cSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE);
255a7a1495cSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE);
256ef66eb69SBarry Smith   } else {
257ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259ef66eb69SBarry Smith     if (*A != *B) {
260ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
262ef66eb69SBarry Smith     }
263ef66eb69SBarry Smith   }
264a7a1495cSBarry Smith   PetscFunctionReturn(0);
265a7a1495cSBarry Smith }
266a7a1495cSBarry Smith 
2674a2ae208SSatish Balay #undef __FUNCT__
2684a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
269d763cef2SBarry Smith /*
270d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
271d763cef2SBarry Smith 
272d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
273d763cef2SBarry Smith    this routine applies the matrix.
274d763cef2SBarry Smith */
27587828ca2SBarry Smith int TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
276d763cef2SBarry Smith {
277d763cef2SBarry Smith   int ierr;
278d763cef2SBarry Smith 
279d763cef2SBarry Smith   PetscFunctionBegin;
280d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
2817c922b88SBarry Smith   PetscValidHeader(x);
2827c922b88SBarry Smith   PetscValidHeader(y);
283d763cef2SBarry Smith 
284d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
285000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
286d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
287000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
288d763cef2SBarry Smith     PetscStackPop;
2897c922b88SBarry Smith   } else {
290000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
291d763cef2SBarry Smith       MatStructure flg;
292d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
293000e7ae3SMatthew Knepley       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
294d763cef2SBarry Smith       PetscStackPop;
295d763cef2SBarry Smith     }
296d763cef2SBarry Smith     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
2977c922b88SBarry Smith   }
298d763cef2SBarry Smith 
299d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
300d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
301d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
302d763cef2SBarry Smith 
303d763cef2SBarry Smith   PetscFunctionReturn(0);
304d763cef2SBarry Smith }
305d763cef2SBarry Smith 
3064a2ae208SSatish Balay #undef __FUNCT__
3074a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
308d763cef2SBarry Smith /*@C
309d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
310d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
311d763cef2SBarry Smith 
312d763cef2SBarry Smith     Collective on TS
313d763cef2SBarry Smith 
314d763cef2SBarry Smith     Input Parameters:
315d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
316d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
317d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
318d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
319d763cef2SBarry Smith 
320d763cef2SBarry Smith     Calling sequence of func:
32187828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
322d763cef2SBarry Smith 
323d763cef2SBarry Smith +   t - current timestep
324d763cef2SBarry Smith .   u - input vector
325d763cef2SBarry Smith .   F - function vector
326d763cef2SBarry Smith -   ctx - [optional] user-defined function context
327d763cef2SBarry Smith 
328d763cef2SBarry Smith     Important:
329d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
330d763cef2SBarry Smith 
331d763cef2SBarry Smith     Level: beginner
332d763cef2SBarry Smith 
333d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
334d763cef2SBarry Smith 
335d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
336d763cef2SBarry Smith @*/
33787828ca2SBarry Smith int TSSetRHSFunction(TS ts,int (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
338d763cef2SBarry Smith {
339d763cef2SBarry Smith   PetscFunctionBegin;
340d763cef2SBarry Smith 
341d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
342d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
34329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
344d763cef2SBarry Smith   }
345000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
346d763cef2SBarry Smith   ts->funP             = ctx;
347d763cef2SBarry Smith   PetscFunctionReturn(0);
348d763cef2SBarry Smith }
349d763cef2SBarry Smith 
3504a2ae208SSatish Balay #undef __FUNCT__
3514a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
352d763cef2SBarry Smith /*@C
353d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
354d763cef2SBarry Smith    Also sets the location to store A.
355d763cef2SBarry Smith 
356d763cef2SBarry Smith    Collective on TS
357d763cef2SBarry Smith 
358d763cef2SBarry Smith    Input Parameters:
359d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
360d763cef2SBarry Smith .  A   - matrix
361d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
362d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
363d763cef2SBarry Smith          if A is not a function of t.
364d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
365d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
366d763cef2SBarry Smith 
367d763cef2SBarry Smith    Calling sequence of func:
36887828ca2SBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,int *flag,void *ctx);
369d763cef2SBarry Smith 
370d763cef2SBarry Smith +  t - current timestep
371d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
372d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
373d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
374d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
375d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
376d763cef2SBarry Smith 
377d763cef2SBarry Smith    Notes:
378d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
379d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
380d763cef2SBarry Smith 
381d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
382d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
383d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
384d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
385d763cef2SBarry Smith    throughout the global iterations.
386d763cef2SBarry Smith 
387d763cef2SBarry Smith    Important:
388d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
389d763cef2SBarry Smith 
390d763cef2SBarry Smith    Level: beginner
391d763cef2SBarry Smith 
392d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
393d763cef2SBarry Smith 
394d763cef2SBarry Smith .seealso: TSSetRHSFunction()
395d763cef2SBarry Smith @*/
39687828ca2SBarry Smith int TSSetRHSMatrix(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
397d763cef2SBarry Smith {
398d763cef2SBarry Smith   PetscFunctionBegin;
399d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
400184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
401184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
402184914b5SBarry Smith   PetscCheckSameComm(ts,A);
403184914b5SBarry Smith   PetscCheckSameComm(ts,B);
404d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
40529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
406d763cef2SBarry Smith   }
407d763cef2SBarry Smith 
408000e7ae3SMatthew Knepley   ts->ops->rhsmatrix = f;
409d763cef2SBarry Smith   ts->jacP           = ctx;
410d763cef2SBarry Smith   ts->A              = A;
411d763cef2SBarry Smith   ts->B              = B;
412d763cef2SBarry Smith 
413d763cef2SBarry Smith   PetscFunctionReturn(0);
414d763cef2SBarry Smith }
415d763cef2SBarry Smith 
4164a2ae208SSatish Balay #undef __FUNCT__
4174a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
418d763cef2SBarry Smith /*@C
419d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
420d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
421d763cef2SBarry Smith 
422d763cef2SBarry Smith    Collective on TS
423d763cef2SBarry Smith 
424d763cef2SBarry Smith    Input Parameters:
425d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
426d763cef2SBarry Smith .  A   - Jacobian matrix
427d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
428d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
429d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
430d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
431d763cef2SBarry Smith 
432d763cef2SBarry Smith    Calling sequence of func:
43387828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
434d763cef2SBarry Smith 
435d763cef2SBarry Smith +  t - current timestep
436d763cef2SBarry Smith .  u - input vector
437d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
438d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
439d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
440d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
441d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
442d763cef2SBarry Smith 
443d763cef2SBarry Smith    Notes:
444d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
445d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
446d763cef2SBarry Smith 
447d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
448d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
449d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
450d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
451d763cef2SBarry Smith    throughout the global iterations.
452d763cef2SBarry Smith 
453d763cef2SBarry Smith    Level: beginner
454d763cef2SBarry Smith 
455d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
456d763cef2SBarry Smith 
457d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
458d763cef2SBarry Smith           SNESDefaultComputeJacobianColor()
459d763cef2SBarry Smith 
460d763cef2SBarry Smith @*/
46187828ca2SBarry Smith int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
462d763cef2SBarry Smith {
463d763cef2SBarry Smith   PetscFunctionBegin;
464d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
465184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
466184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
467184914b5SBarry Smith   PetscCheckSameComm(ts,A);
468184914b5SBarry Smith   PetscCheckSameComm(ts,B);
469d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
47029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
471d763cef2SBarry Smith   }
472d763cef2SBarry Smith 
473000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
474d763cef2SBarry Smith   ts->jacP             = ctx;
475d763cef2SBarry Smith   ts->A                = A;
476d763cef2SBarry Smith   ts->B                = B;
477d763cef2SBarry Smith   PetscFunctionReturn(0);
478d763cef2SBarry Smith }
479d763cef2SBarry Smith 
4804a2ae208SSatish Balay #undef __FUNCT__
4814a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
482d763cef2SBarry Smith /*
483d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
484d763cef2SBarry Smith 
485d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
486d763cef2SBarry Smith    this routine applies the matrix.
487d763cef2SBarry Smith */
48887828ca2SBarry Smith int TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
489d763cef2SBarry Smith {
490d763cef2SBarry Smith   int ierr;
491d763cef2SBarry Smith 
492d763cef2SBarry Smith   PetscFunctionBegin;
493d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
494d763cef2SBarry Smith   PetscValidHeader(x);
495184914b5SBarry Smith   PetscCheckSameComm(ts,x);
496d763cef2SBarry Smith 
497000e7ae3SMatthew Knepley   if (ts->ops->rhsbc) {
498d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
499000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
500d763cef2SBarry Smith     PetscStackPop;
501d763cef2SBarry Smith     PetscFunctionReturn(0);
502d763cef2SBarry Smith   }
503d763cef2SBarry Smith 
504d763cef2SBarry Smith   PetscFunctionReturn(0);
505d763cef2SBarry Smith }
506d763cef2SBarry Smith 
5074a2ae208SSatish Balay #undef __FUNCT__
5084a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
509d763cef2SBarry Smith /*@C
510d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
511d763cef2SBarry Smith     boundary conditions for the function F.
512d763cef2SBarry Smith 
513d763cef2SBarry Smith     Collective on TS
514d763cef2SBarry Smith 
515d763cef2SBarry Smith     Input Parameters:
516d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
517d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
518d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
519d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
520d763cef2SBarry Smith 
521d763cef2SBarry Smith     Calling sequence of func:
52287828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
523d763cef2SBarry Smith 
524d763cef2SBarry Smith +   t - current timestep
525d763cef2SBarry Smith .   F - function vector
526d763cef2SBarry Smith -   ctx - [optional] user-defined function context
527d763cef2SBarry Smith 
528d763cef2SBarry Smith     Level: intermediate
529d763cef2SBarry Smith 
530d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
531d763cef2SBarry Smith @*/
53287828ca2SBarry Smith int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,PetscReal,Vec,void*),void *ctx)
533d763cef2SBarry Smith {
534d763cef2SBarry Smith   PetscFunctionBegin;
535d763cef2SBarry Smith 
536d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
537d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
53829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
539d763cef2SBarry Smith   }
540000e7ae3SMatthew Knepley   ts->ops->rhsbc = f;
541d763cef2SBarry Smith   ts->bcP        = ctx;
542d763cef2SBarry Smith   PetscFunctionReturn(0);
543d763cef2SBarry Smith }
544d763cef2SBarry Smith 
5454a2ae208SSatish Balay #undef __FUNCT__
5464a2ae208SSatish Balay #define __FUNCT__ "TSView"
5477e2c5f70SBarry Smith /*@C
548d763cef2SBarry Smith     TSView - Prints the TS data structure.
549d763cef2SBarry Smith 
5504c49b128SBarry Smith     Collective on TS
551d763cef2SBarry Smith 
552d763cef2SBarry Smith     Input Parameters:
553d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
554d763cef2SBarry Smith -   viewer - visualization context
555d763cef2SBarry Smith 
556d763cef2SBarry Smith     Options Database Key:
557d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
558d763cef2SBarry Smith 
559d763cef2SBarry Smith     Notes:
560d763cef2SBarry Smith     The available visualization contexts include
561b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
562b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
563d763cef2SBarry Smith          output where only the first processor opens
564d763cef2SBarry Smith          the file.  All other processors send their
565d763cef2SBarry Smith          data to the first processor to print.
566d763cef2SBarry Smith 
567d763cef2SBarry Smith     The user can open an alternative visualization context with
568b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
569d763cef2SBarry Smith 
570d763cef2SBarry Smith     Level: beginner
571d763cef2SBarry Smith 
572d763cef2SBarry Smith .keywords: TS, timestep, view
573d763cef2SBarry Smith 
574b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
575d763cef2SBarry Smith @*/
576b0a32e0cSBarry Smith int TSView(TS ts,PetscViewer viewer)
577d763cef2SBarry Smith {
578d763cef2SBarry Smith   int        ierr;
579454a90a3SBarry Smith   char       *type;
5806831982aSBarry Smith   PetscTruth isascii,isstring;
581d763cef2SBarry Smith 
582d763cef2SBarry Smith   PetscFunctionBegin;
583d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
584b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
585b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
5866831982aSBarry Smith   PetscCheckSameComm(ts,viewer);
587fd16b177SBarry Smith 
588b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
589b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
5900f5bd95cSBarry Smith   if (isascii) {
591b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
592454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
593454a90a3SBarry Smith     if (type) {
594b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
595184914b5SBarry Smith     } else {
596b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
597184914b5SBarry Smith     }
598000e7ae3SMatthew Knepley     if (ts->ops->view) {
599b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
600000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
601b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
602d763cef2SBarry Smith     }
603b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
604b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
605d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
606b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
607d763cef2SBarry Smith     }
608b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
6090f5bd95cSBarry Smith   } else if (isstring) {
610454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
611b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
612d763cef2SBarry Smith   }
613b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
614d763cef2SBarry Smith   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
615d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
616b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
617d763cef2SBarry Smith   PetscFunctionReturn(0);
618d763cef2SBarry Smith }
619d763cef2SBarry Smith 
620d763cef2SBarry Smith 
6214a2ae208SSatish Balay #undef __FUNCT__
6224a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
623d763cef2SBarry Smith /*@C
624d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
625d763cef2SBarry Smith    the timesteppers.
626d763cef2SBarry Smith 
627d763cef2SBarry Smith    Collective on TS
628d763cef2SBarry Smith 
629d763cef2SBarry Smith    Input Parameters:
630d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
631d763cef2SBarry Smith -  usrP - optional user context
632d763cef2SBarry Smith 
633d763cef2SBarry Smith    Level: intermediate
634d763cef2SBarry Smith 
635d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
636d763cef2SBarry Smith 
637d763cef2SBarry Smith .seealso: TSGetApplicationContext()
638d763cef2SBarry Smith @*/
639d763cef2SBarry Smith int TSSetApplicationContext(TS ts,void *usrP)
640d763cef2SBarry Smith {
641d763cef2SBarry Smith   PetscFunctionBegin;
642d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
643d763cef2SBarry Smith   ts->user = usrP;
644d763cef2SBarry Smith   PetscFunctionReturn(0);
645d763cef2SBarry Smith }
646d763cef2SBarry Smith 
6474a2ae208SSatish Balay #undef __FUNCT__
6484a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
649d763cef2SBarry Smith /*@C
650d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
651d763cef2SBarry Smith     timestepper.
652d763cef2SBarry Smith 
653d763cef2SBarry Smith     Not Collective
654d763cef2SBarry Smith 
655d763cef2SBarry Smith     Input Parameter:
656d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
657d763cef2SBarry Smith 
658d763cef2SBarry Smith     Output Parameter:
659d763cef2SBarry Smith .   usrP - user context
660d763cef2SBarry Smith 
661d763cef2SBarry Smith     Level: intermediate
662d763cef2SBarry Smith 
663d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
664d763cef2SBarry Smith 
665d763cef2SBarry Smith .seealso: TSSetApplicationContext()
666d763cef2SBarry Smith @*/
667d763cef2SBarry Smith int TSGetApplicationContext(TS ts,void **usrP)
668d763cef2SBarry Smith {
669d763cef2SBarry Smith   PetscFunctionBegin;
670d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
671d763cef2SBarry Smith   *usrP = ts->user;
672d763cef2SBarry Smith   PetscFunctionReturn(0);
673d763cef2SBarry Smith }
674d763cef2SBarry Smith 
6754a2ae208SSatish Balay #undef __FUNCT__
6764a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
677d763cef2SBarry Smith /*@
678d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
679d763cef2SBarry Smith 
680d763cef2SBarry Smith    Not Collective
681d763cef2SBarry Smith 
682d763cef2SBarry Smith    Input Parameter:
683d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
684d763cef2SBarry Smith 
685d763cef2SBarry Smith    Output Parameter:
686d763cef2SBarry Smith .  iter - number steps so far
687d763cef2SBarry Smith 
688d763cef2SBarry Smith    Level: intermediate
689d763cef2SBarry Smith 
690d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
691d763cef2SBarry Smith @*/
692d763cef2SBarry Smith int TSGetTimeStepNumber(TS ts,int* iter)
693d763cef2SBarry Smith {
694d763cef2SBarry Smith   PetscFunctionBegin;
695d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
696d763cef2SBarry Smith   *iter = ts->steps;
697d763cef2SBarry Smith   PetscFunctionReturn(0);
698d763cef2SBarry Smith }
699d763cef2SBarry Smith 
7004a2ae208SSatish Balay #undef __FUNCT__
7014a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
702d763cef2SBarry Smith /*@
703d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
704d763cef2SBarry Smith    as well as the initial time.
705d763cef2SBarry Smith 
706d763cef2SBarry Smith    Collective on TS
707d763cef2SBarry Smith 
708d763cef2SBarry Smith    Input Parameters:
709d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
710d763cef2SBarry Smith .  initial_time - the initial time
711d763cef2SBarry Smith -  time_step - the size of the timestep
712d763cef2SBarry Smith 
713d763cef2SBarry Smith    Level: intermediate
714d763cef2SBarry Smith 
715d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
716d763cef2SBarry Smith 
717d763cef2SBarry Smith .keywords: TS, set, initial, timestep
718d763cef2SBarry Smith @*/
71987828ca2SBarry Smith int TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
720d763cef2SBarry Smith {
721d763cef2SBarry Smith   PetscFunctionBegin;
722d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
723d763cef2SBarry Smith   ts->time_step         = time_step;
724d763cef2SBarry Smith   ts->initial_time_step = time_step;
725d763cef2SBarry Smith   ts->ptime             = initial_time;
726d763cef2SBarry Smith   PetscFunctionReturn(0);
727d763cef2SBarry Smith }
728d763cef2SBarry Smith 
7294a2ae208SSatish Balay #undef __FUNCT__
7304a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
731d763cef2SBarry Smith /*@
732d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
733d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
734d763cef2SBarry Smith 
735d763cef2SBarry Smith    Collective on TS
736d763cef2SBarry Smith 
737d763cef2SBarry Smith    Input Parameters:
738d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
739d763cef2SBarry Smith -  time_step - the size of the timestep
740d763cef2SBarry Smith 
741d763cef2SBarry Smith    Level: intermediate
742d763cef2SBarry Smith 
743d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
744d763cef2SBarry Smith 
745d763cef2SBarry Smith .keywords: TS, set, timestep
746d763cef2SBarry Smith @*/
74787828ca2SBarry Smith int TSSetTimeStep(TS ts,PetscReal time_step)
748d763cef2SBarry Smith {
749d763cef2SBarry Smith   PetscFunctionBegin;
750d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
751d763cef2SBarry Smith   ts->time_step = time_step;
752d763cef2SBarry Smith   PetscFunctionReturn(0);
753d763cef2SBarry Smith }
754d763cef2SBarry Smith 
7554a2ae208SSatish Balay #undef __FUNCT__
7564a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
757d763cef2SBarry Smith /*@
758d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
759d763cef2SBarry Smith 
760d763cef2SBarry Smith    Not Collective
761d763cef2SBarry Smith 
762d763cef2SBarry Smith    Input Parameter:
763d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
764d763cef2SBarry Smith 
765d763cef2SBarry Smith    Output Parameter:
766d763cef2SBarry Smith .  dt - the current timestep size
767d763cef2SBarry Smith 
768d763cef2SBarry Smith    Level: intermediate
769d763cef2SBarry Smith 
770d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
771d763cef2SBarry Smith 
772d763cef2SBarry Smith .keywords: TS, get, timestep
773d763cef2SBarry Smith @*/
77487828ca2SBarry Smith int TSGetTimeStep(TS ts,PetscReal* dt)
775d763cef2SBarry Smith {
776d763cef2SBarry Smith   PetscFunctionBegin;
777d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
778d763cef2SBarry Smith   *dt = ts->time_step;
779d763cef2SBarry Smith   PetscFunctionReturn(0);
780d763cef2SBarry Smith }
781d763cef2SBarry Smith 
7824a2ae208SSatish Balay #undef __FUNCT__
7834a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
784d763cef2SBarry Smith /*@C
785d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
786d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
787d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
788d763cef2SBarry Smith    the solution at the next timestep has been calculated.
789d763cef2SBarry Smith 
790d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
791d763cef2SBarry Smith 
792d763cef2SBarry Smith    Input Parameter:
793d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
794d763cef2SBarry Smith 
795d763cef2SBarry Smith    Output Parameter:
796d763cef2SBarry Smith .  v - the vector containing the solution
797d763cef2SBarry Smith 
798d763cef2SBarry Smith    Level: intermediate
799d763cef2SBarry Smith 
800d763cef2SBarry Smith .seealso: TSGetTimeStep()
801d763cef2SBarry Smith 
802d763cef2SBarry Smith .keywords: TS, timestep, get, solution
803d763cef2SBarry Smith @*/
804d763cef2SBarry Smith int TSGetSolution(TS ts,Vec *v)
805d763cef2SBarry Smith {
806d763cef2SBarry Smith   PetscFunctionBegin;
807d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
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;
836bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
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;
866bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
867bdad233fSMatthew Knepley   PetscValidPointer(type);
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;
901d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
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;
933d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
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 
939d763cef2SBarry Smith   if (ts->sles) {ierr = SLESDestroy(ts->sles);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
969d763cef2SBarry Smith    SLES, 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;
981d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
98229bbc08cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetSLES()");
983d763cef2SBarry Smith   *snes = ts->snes;
984d763cef2SBarry Smith   PetscFunctionReturn(0);
985d763cef2SBarry Smith }
986d763cef2SBarry Smith 
9874a2ae208SSatish Balay #undef __FUNCT__
9884a2ae208SSatish Balay #define __FUNCT__ "TSGetSLES"
989d763cef2SBarry Smith /*@C
990d763cef2SBarry Smith    TSGetSLES - Returns the SLES (linear solver) associated with
991d763cef2SBarry Smith    a TS (timestepper) context.
992d763cef2SBarry Smith 
993d763cef2SBarry Smith    Not Collective, but SLES is parallel if TS is parallel
994d763cef2SBarry Smith 
995d763cef2SBarry Smith    Input Parameter:
996d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
997d763cef2SBarry Smith 
998d763cef2SBarry Smith    Output Parameter:
999d763cef2SBarry Smith .  sles - the nonlinear solver context
1000d763cef2SBarry Smith 
1001d763cef2SBarry Smith    Notes:
1002d763cef2SBarry Smith    The user can then directly manipulate the SLES context to set various
1003d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1004d763cef2SBarry Smith    KSP and PC contexts as well.
1005d763cef2SBarry Smith 
1006d763cef2SBarry Smith    TSGetSLES() does not work for integrators that do not use SLES;
1007d763cef2SBarry Smith    in this case TSGetSLES() returns PETSC_NULL in sles.
1008d763cef2SBarry Smith 
1009d763cef2SBarry Smith    Level: beginner
1010d763cef2SBarry Smith 
1011d763cef2SBarry Smith .keywords: timestep, get, SLES
1012d763cef2SBarry Smith @*/
1013d763cef2SBarry Smith int TSGetSLES(TS ts,SLES *sles)
1014d763cef2SBarry Smith {
1015d763cef2SBarry Smith   PetscFunctionBegin;
1016d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
101729bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1018d763cef2SBarry Smith   *sles = ts->sles;
1019d763cef2SBarry Smith   PetscFunctionReturn(0);
1020d763cef2SBarry Smith }
1021d763cef2SBarry Smith 
1022d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1023d763cef2SBarry Smith 
10244a2ae208SSatish Balay #undef __FUNCT__
10254a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1026d763cef2SBarry Smith /*@
1027d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1028d763cef2SBarry Smith    maximum time for iteration.
1029d763cef2SBarry Smith 
1030d763cef2SBarry Smith    Collective on TS
1031d763cef2SBarry Smith 
1032d763cef2SBarry Smith    Input Parameters:
1033d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1034d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1035d763cef2SBarry Smith -  maxtime - final time to iterate to
1036d763cef2SBarry Smith 
1037d763cef2SBarry Smith    Options Database Keys:
1038d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1039d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1040d763cef2SBarry Smith 
1041d763cef2SBarry Smith    Notes:
1042d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1043d763cef2SBarry Smith 
1044d763cef2SBarry Smith    Level: intermediate
1045d763cef2SBarry Smith 
1046d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1047d763cef2SBarry Smith @*/
104887828ca2SBarry Smith int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime)
1049d763cef2SBarry Smith {
1050d763cef2SBarry Smith   PetscFunctionBegin;
1051d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1052d763cef2SBarry Smith   ts->max_steps = maxsteps;
1053d763cef2SBarry Smith   ts->max_time  = maxtime;
1054d763cef2SBarry Smith   PetscFunctionReturn(0);
1055d763cef2SBarry Smith }
1056d763cef2SBarry Smith 
10574a2ae208SSatish Balay #undef __FUNCT__
10584a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1059d763cef2SBarry Smith /*@
1060d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1061d763cef2SBarry Smith    for use by the TS routines.
1062d763cef2SBarry Smith 
1063d763cef2SBarry Smith    Collective on TS and Vec
1064d763cef2SBarry Smith 
1065d763cef2SBarry Smith    Input Parameters:
1066d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1067d763cef2SBarry Smith -  x - the solution vector
1068d763cef2SBarry Smith 
1069d763cef2SBarry Smith    Level: beginner
1070d763cef2SBarry Smith 
1071d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1072d763cef2SBarry Smith @*/
1073d763cef2SBarry Smith int TSSetSolution(TS ts,Vec x)
1074d763cef2SBarry Smith {
1075d763cef2SBarry Smith   PetscFunctionBegin;
1076d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1077d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1078d763cef2SBarry Smith   PetscFunctionReturn(0);
1079d763cef2SBarry Smith }
1080d763cef2SBarry Smith 
1081e74ef692SMatthew Knepley #undef __FUNCT__
1082e74ef692SMatthew Knepley #define __FUNCT__ "TSSetRhsBC"
1083000e7ae3SMatthew Knepley /*@
1084000e7ae3SMatthew Knepley   TSSetRhsBC - Sets the function which applies boundary conditions
1085000e7ae3SMatthew Knepley   to the Rhs of each system.
1086000e7ae3SMatthew Knepley 
1087000e7ae3SMatthew Knepley   Collective on TS
1088000e7ae3SMatthew Knepley 
1089000e7ae3SMatthew Knepley   Input Parameters:
1090000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1091000e7ae3SMatthew Knepley - func - The function
1092000e7ae3SMatthew Knepley 
1093000e7ae3SMatthew Knepley   Calling sequence of func:
1094000e7ae3SMatthew Knepley . func (TS ts, Vec rhs, void *ctx);
1095000e7ae3SMatthew Knepley 
1096000e7ae3SMatthew Knepley + rhs - The current rhs vector
1097000e7ae3SMatthew Knepley - ctx - The user-context
1098000e7ae3SMatthew Knepley 
1099000e7ae3SMatthew Knepley   Level: intermediate
1100000e7ae3SMatthew Knepley 
1101000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1102000e7ae3SMatthew Knepley @*/
1103000e7ae3SMatthew Knepley int TSSetRhsBC(TS ts, int (*func)(TS, Vec, void *))
1104000e7ae3SMatthew Knepley {
1105000e7ae3SMatthew Knepley   PetscFunctionBegin;
1106000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1107000e7ae3SMatthew Knepley   ts->ops->applyrhsbc = func;
1108000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1109000e7ae3SMatthew Knepley }
1110000e7ae3SMatthew Knepley 
1111e74ef692SMatthew Knepley #undef __FUNCT__
1112e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultRhsBC"
1113000e7ae3SMatthew Knepley /*@
1114000e7ae3SMatthew Knepley   TSDefaultRhsBC - The default boundary condition function which does nothing.
1115000e7ae3SMatthew Knepley 
1116000e7ae3SMatthew Knepley   Collective on TS
1117000e7ae3SMatthew Knepley 
1118000e7ae3SMatthew Knepley   Input Parameters:
1119000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1120000e7ae3SMatthew Knepley . rhs - The Rhs
1121000e7ae3SMatthew Knepley - ctx - The user-context
1122000e7ae3SMatthew Knepley 
1123000e7ae3SMatthew Knepley   Level: developer
1124000e7ae3SMatthew Knepley 
1125000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1126000e7ae3SMatthew Knepley @*/
1127000e7ae3SMatthew Knepley int TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1128000e7ae3SMatthew Knepley {
1129000e7ae3SMatthew Knepley   PetscFunctionBegin;
1130000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1131000e7ae3SMatthew Knepley }
1132000e7ae3SMatthew Knepley 
1133e74ef692SMatthew Knepley #undef __FUNCT__
1134e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSystemMatrixBC"
1135000e7ae3SMatthew Knepley /*@
1136000e7ae3SMatthew Knepley   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1137000e7ae3SMatthew Knepley   to the system matrix and preconditioner of each system.
1138000e7ae3SMatthew Knepley 
1139000e7ae3SMatthew Knepley   Collective on TS
1140000e7ae3SMatthew Knepley 
1141000e7ae3SMatthew Knepley   Input Parameters:
1142000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1143000e7ae3SMatthew Knepley - func - The function
1144000e7ae3SMatthew Knepley 
1145000e7ae3SMatthew Knepley   Calling sequence of func:
1146000e7ae3SMatthew Knepley . func (TS ts, Mat A, Mat B, void *ctx);
1147000e7ae3SMatthew Knepley 
1148000e7ae3SMatthew Knepley + A   - The current system matrix
1149000e7ae3SMatthew Knepley . B   - The current preconditioner
1150000e7ae3SMatthew Knepley - ctx - The user-context
1151000e7ae3SMatthew Knepley 
1152000e7ae3SMatthew Knepley   Level: intermediate
1153000e7ae3SMatthew Knepley 
1154000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1155000e7ae3SMatthew Knepley @*/
1156000e7ae3SMatthew Knepley int TSSetSystemMatrixBC(TS ts, int (*func)(TS, Mat, Mat, void *))
1157000e7ae3SMatthew Knepley {
1158000e7ae3SMatthew Knepley   PetscFunctionBegin;
1159000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1160000e7ae3SMatthew Knepley   ts->ops->applymatrixbc = func;
1161000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1162000e7ae3SMatthew Knepley }
1163000e7ae3SMatthew Knepley 
1164e74ef692SMatthew Knepley #undef __FUNCT__
1165e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSystemMatrixBC"
1166000e7ae3SMatthew Knepley /*@
1167000e7ae3SMatthew Knepley   TSDefaultSystemMatrixBC - The default boundary condition function which
1168000e7ae3SMatthew Knepley   does nothing.
1169000e7ae3SMatthew Knepley 
1170000e7ae3SMatthew Knepley   Collective on TS
1171000e7ae3SMatthew Knepley 
1172000e7ae3SMatthew Knepley   Input Parameters:
1173000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1174000e7ae3SMatthew Knepley . A   - The system matrix
1175000e7ae3SMatthew Knepley . B   - The preconditioner
1176000e7ae3SMatthew Knepley - ctx - The user-context
1177000e7ae3SMatthew Knepley 
1178000e7ae3SMatthew Knepley   Level: developer
1179000e7ae3SMatthew Knepley 
1180000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1181000e7ae3SMatthew Knepley @*/
1182000e7ae3SMatthew Knepley int TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1183000e7ae3SMatthew Knepley {
1184000e7ae3SMatthew Knepley   PetscFunctionBegin;
1185000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1186000e7ae3SMatthew Knepley }
1187000e7ae3SMatthew Knepley 
1188e74ef692SMatthew Knepley #undef __FUNCT__
1189e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSolutionBC"
1190000e7ae3SMatthew Knepley /*@
1191000e7ae3SMatthew Knepley   TSSetSolutionBC - Sets the function which applies boundary conditions
1192000e7ae3SMatthew Knepley   to the solution of each system. This is necessary in nonlinear systems
1193000e7ae3SMatthew Knepley   which time dependent boundary conditions.
1194000e7ae3SMatthew Knepley 
1195000e7ae3SMatthew Knepley   Collective on TS
1196000e7ae3SMatthew Knepley 
1197000e7ae3SMatthew Knepley   Input Parameters:
1198000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1199000e7ae3SMatthew Knepley - func - The function
1200000e7ae3SMatthew Knepley 
1201000e7ae3SMatthew Knepley   Calling sequence of func:
1202000e7ae3SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
1203000e7ae3SMatthew Knepley 
1204000e7ae3SMatthew Knepley + sol - The current solution vector
1205000e7ae3SMatthew Knepley - ctx - The user-context
1206000e7ae3SMatthew Knepley 
1207000e7ae3SMatthew Knepley   Level: intermediate
1208000e7ae3SMatthew Knepley 
1209000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1210000e7ae3SMatthew Knepley @*/
1211000e7ae3SMatthew Knepley int TSSetSolutionBC(TS ts, int (*func)(TS, Vec, void *))
1212000e7ae3SMatthew Knepley {
1213000e7ae3SMatthew Knepley   PetscFunctionBegin;
1214000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1215000e7ae3SMatthew Knepley   ts->ops->applysolbc = func;
1216000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1217000e7ae3SMatthew Knepley }
1218000e7ae3SMatthew Knepley 
1219e74ef692SMatthew Knepley #undef __FUNCT__
1220e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSolutionBC"
1221000e7ae3SMatthew Knepley /*@
1222000e7ae3SMatthew Knepley   TSDefaultSolutionBC - The default boundary condition function which
1223000e7ae3SMatthew Knepley   does nothing.
1224000e7ae3SMatthew Knepley 
1225000e7ae3SMatthew Knepley   Collective on TS
1226000e7ae3SMatthew Knepley 
1227000e7ae3SMatthew Knepley   Input Parameters:
1228000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1229000e7ae3SMatthew Knepley . sol - The solution
1230000e7ae3SMatthew Knepley - ctx - The user-context
1231000e7ae3SMatthew Knepley 
1232000e7ae3SMatthew Knepley   Level: developer
1233000e7ae3SMatthew Knepley 
1234000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1235000e7ae3SMatthew Knepley @*/
1236000e7ae3SMatthew Knepley int TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1237000e7ae3SMatthew Knepley {
1238000e7ae3SMatthew Knepley   PetscFunctionBegin;
1239000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1240000e7ae3SMatthew Knepley }
1241000e7ae3SMatthew Knepley 
1242e74ef692SMatthew Knepley #undef __FUNCT__
1243e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1244000e7ae3SMatthew Knepley /*@
1245000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1246000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1247000e7ae3SMatthew Knepley 
1248000e7ae3SMatthew Knepley   Collective on TS
1249000e7ae3SMatthew Knepley 
1250000e7ae3SMatthew Knepley   Input Parameters:
1251000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1252000e7ae3SMatthew Knepley - func - The function
1253000e7ae3SMatthew Knepley 
1254000e7ae3SMatthew Knepley   Calling sequence of func:
1255000e7ae3SMatthew Knepley . func (TS ts);
1256000e7ae3SMatthew Knepley 
1257000e7ae3SMatthew Knepley   Level: intermediate
1258000e7ae3SMatthew Knepley 
1259000e7ae3SMatthew Knepley .keywords: TS, timestep
1260000e7ae3SMatthew Knepley @*/
1261000e7ae3SMatthew Knepley int TSSetPreStep(TS ts, int (*func)(TS))
1262000e7ae3SMatthew Knepley {
1263000e7ae3SMatthew Knepley   PetscFunctionBegin;
1264000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1265000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1266000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1267000e7ae3SMatthew Knepley }
1268000e7ae3SMatthew Knepley 
1269e74ef692SMatthew Knepley #undef __FUNCT__
1270e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1271000e7ae3SMatthew Knepley /*@
1272000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1273000e7ae3SMatthew Knepley 
1274000e7ae3SMatthew Knepley   Collective on TS
1275000e7ae3SMatthew Knepley 
1276000e7ae3SMatthew Knepley   Input Parameters:
1277000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1278000e7ae3SMatthew Knepley 
1279000e7ae3SMatthew Knepley   Level: developer
1280000e7ae3SMatthew Knepley 
1281000e7ae3SMatthew Knepley .keywords: TS, timestep
1282000e7ae3SMatthew Knepley @*/
1283000e7ae3SMatthew Knepley int TSDefaultPreStep(TS ts)
1284000e7ae3SMatthew Knepley {
1285000e7ae3SMatthew Knepley   PetscFunctionBegin;
1286000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1287000e7ae3SMatthew Knepley }
1288000e7ae3SMatthew Knepley 
1289e74ef692SMatthew Knepley #undef __FUNCT__
1290e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1291000e7ae3SMatthew Knepley /*@
1292000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1293000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1294000e7ae3SMatthew Knepley   the time step.
1295000e7ae3SMatthew Knepley 
1296000e7ae3SMatthew Knepley   Collective on TS
1297000e7ae3SMatthew Knepley 
1298000e7ae3SMatthew Knepley   Input Parameters:
1299000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1300000e7ae3SMatthew Knepley - func - The function
1301000e7ae3SMatthew Knepley 
1302000e7ae3SMatthew Knepley   Calling sequence of func:
1303000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1304000e7ae3SMatthew Knepley 
1305000e7ae3SMatthew Knepley + t   - The current time
1306000e7ae3SMatthew Knepley - dt  - The current time step
1307000e7ae3SMatthew Knepley 
1308000e7ae3SMatthew Knepley   Level: intermediate
1309000e7ae3SMatthew Knepley 
1310000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1311000e7ae3SMatthew Knepley @*/
1312000e7ae3SMatthew Knepley int TSSetUpdate(TS ts, int (*func)(TS, double, double *))
1313000e7ae3SMatthew Knepley {
1314000e7ae3SMatthew Knepley   PetscFunctionBegin;
1315000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1316000e7ae3SMatthew Knepley   ts->ops->update = func;
1317000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1318000e7ae3SMatthew Knepley }
1319000e7ae3SMatthew Knepley 
1320e74ef692SMatthew Knepley #undef __FUNCT__
1321e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1322000e7ae3SMatthew Knepley /*@
1323000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1324000e7ae3SMatthew Knepley 
1325000e7ae3SMatthew Knepley   Collective on TS
1326000e7ae3SMatthew Knepley 
1327000e7ae3SMatthew Knepley   Input Parameters:
1328000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1329000e7ae3SMatthew Knepley - t   - The current time
1330000e7ae3SMatthew Knepley 
1331000e7ae3SMatthew Knepley   Output Parameters:
1332000e7ae3SMatthew Knepley . dt  - The current time step
1333000e7ae3SMatthew Knepley 
1334000e7ae3SMatthew Knepley   Level: developer
1335000e7ae3SMatthew Knepley 
1336000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1337000e7ae3SMatthew Knepley @*/
1338000e7ae3SMatthew Knepley int TSDefaultUpdate(TS ts, double t, double *dt)
1339000e7ae3SMatthew Knepley {
1340000e7ae3SMatthew Knepley   PetscFunctionBegin;
1341000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1342000e7ae3SMatthew Knepley }
1343000e7ae3SMatthew Knepley 
1344e74ef692SMatthew Knepley #undef __FUNCT__
1345e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1346000e7ae3SMatthew Knepley /*@
1347000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1348000e7ae3SMatthew Knepley   called once at the end of time stepping.
1349000e7ae3SMatthew Knepley 
1350000e7ae3SMatthew Knepley   Collective on TS
1351000e7ae3SMatthew Knepley 
1352000e7ae3SMatthew Knepley   Input Parameters:
1353000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1354000e7ae3SMatthew Knepley - func - The function
1355000e7ae3SMatthew Knepley 
1356000e7ae3SMatthew Knepley   Calling sequence of func:
1357000e7ae3SMatthew Knepley . func (TS ts);
1358000e7ae3SMatthew Knepley 
1359000e7ae3SMatthew Knepley   Level: intermediate
1360000e7ae3SMatthew Knepley 
1361000e7ae3SMatthew Knepley .keywords: TS, timestep
1362000e7ae3SMatthew Knepley @*/
1363000e7ae3SMatthew Knepley int TSSetPostStep(TS ts, int (*func)(TS))
1364000e7ae3SMatthew Knepley {
1365000e7ae3SMatthew Knepley   PetscFunctionBegin;
1366000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1367000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1368000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1369000e7ae3SMatthew Knepley }
1370000e7ae3SMatthew Knepley 
1371e74ef692SMatthew Knepley #undef __FUNCT__
1372e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1373000e7ae3SMatthew Knepley /*@
1374000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1375000e7ae3SMatthew Knepley 
1376000e7ae3SMatthew Knepley   Collective on TS
1377000e7ae3SMatthew Knepley 
1378000e7ae3SMatthew Knepley   Input Parameters:
1379000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1380000e7ae3SMatthew Knepley 
1381000e7ae3SMatthew Knepley   Level: developer
1382000e7ae3SMatthew Knepley 
1383000e7ae3SMatthew Knepley .keywords: TS, timestep
1384000e7ae3SMatthew Knepley @*/
1385000e7ae3SMatthew Knepley int TSDefaultPostStep(TS ts)
1386000e7ae3SMatthew Knepley {
1387000e7ae3SMatthew Knepley   PetscFunctionBegin;
1388000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1389000e7ae3SMatthew Knepley }
1390000e7ae3SMatthew Knepley 
1391d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1392d763cef2SBarry Smith 
13934a2ae208SSatish Balay #undef __FUNCT__
13944a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
1395d763cef2SBarry Smith /*@C
1396d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1397d763cef2SBarry Smith    timestep to display the iteration's  progress.
1398d763cef2SBarry Smith 
1399d763cef2SBarry Smith    Collective on TS
1400d763cef2SBarry Smith 
1401d763cef2SBarry Smith    Input Parameters:
1402d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1403d763cef2SBarry Smith .  func - monitoring routine
1404329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1405b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1406b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1407b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1408d763cef2SBarry Smith 
1409d763cef2SBarry Smith    Calling sequence of func:
141087828ca2SBarry Smith $    int func(TS ts,int steps,PetscReal time,Vec x,void *mctx)
1411d763cef2SBarry Smith 
1412d763cef2SBarry Smith +    ts - the TS context
1413d763cef2SBarry Smith .    steps - iteration number
14141f06c33eSBarry Smith .    time - current time
1415d763cef2SBarry Smith .    x - current iterate
1416d763cef2SBarry Smith -    mctx - [optional] monitoring context
1417d763cef2SBarry Smith 
1418d763cef2SBarry Smith    Notes:
1419d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1420d763cef2SBarry Smith    already has been loaded.
1421d763cef2SBarry Smith 
1422d763cef2SBarry Smith    Level: intermediate
1423d763cef2SBarry Smith 
1424d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1425d763cef2SBarry Smith 
1426d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
1427d763cef2SBarry Smith @*/
142887828ca2SBarry Smith int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*))
1429d763cef2SBarry Smith {
1430d763cef2SBarry Smith   PetscFunctionBegin;
1431d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1432d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
143329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1434d763cef2SBarry Smith   }
1435d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1436329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1437d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1438d763cef2SBarry Smith   PetscFunctionReturn(0);
1439d763cef2SBarry Smith }
1440d763cef2SBarry Smith 
14414a2ae208SSatish Balay #undef __FUNCT__
14424a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
1443d763cef2SBarry Smith /*@C
1444d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1445d763cef2SBarry Smith 
1446d763cef2SBarry Smith    Collective on TS
1447d763cef2SBarry Smith 
1448d763cef2SBarry Smith    Input Parameters:
1449d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1450d763cef2SBarry Smith 
1451d763cef2SBarry Smith    Notes:
1452d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1453d763cef2SBarry Smith 
1454d763cef2SBarry Smith    Level: intermediate
1455d763cef2SBarry Smith 
1456d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1457d763cef2SBarry Smith 
1458d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
1459d763cef2SBarry Smith @*/
1460d763cef2SBarry Smith int TSClearMonitor(TS ts)
1461d763cef2SBarry Smith {
1462d763cef2SBarry Smith   PetscFunctionBegin;
1463d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1464d763cef2SBarry Smith   ts->numbermonitors = 0;
1465d763cef2SBarry Smith   PetscFunctionReturn(0);
1466d763cef2SBarry Smith }
1467d763cef2SBarry Smith 
14684a2ae208SSatish Balay #undef __FUNCT__
14694a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
147087828ca2SBarry Smith int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
1471d763cef2SBarry Smith {
1472d132466eSBarry Smith   int ierr;
1473d132466eSBarry Smith 
1474d763cef2SBarry Smith   PetscFunctionBegin;
1475142b95e3SSatish Balay   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1476d763cef2SBarry Smith   PetscFunctionReturn(0);
1477d763cef2SBarry Smith }
1478d763cef2SBarry Smith 
14794a2ae208SSatish Balay #undef __FUNCT__
14804a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1481d763cef2SBarry Smith /*@
1482d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1483d763cef2SBarry Smith 
1484d763cef2SBarry Smith    Collective on TS
1485d763cef2SBarry Smith 
1486d763cef2SBarry Smith    Input Parameter:
1487d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1488d763cef2SBarry Smith 
1489d763cef2SBarry Smith    Output Parameters:
1490d763cef2SBarry Smith +  steps - number of iterations until termination
1491142b95e3SSatish Balay -  ptime - time until termination
1492d763cef2SBarry Smith 
1493d763cef2SBarry Smith    Level: beginner
1494d763cef2SBarry Smith 
1495d763cef2SBarry Smith .keywords: TS, timestep, solve
1496d763cef2SBarry Smith 
1497d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1498d763cef2SBarry Smith @*/
149987828ca2SBarry Smith int TSStep(TS ts,int *steps,PetscReal *ptime)
1500d763cef2SBarry Smith {
1501d405a339SMatthew Knepley   PetscTruth opt;
1502f1af5d2fSBarry Smith   int        ierr;
1503d763cef2SBarry Smith 
1504d763cef2SBarry Smith   PetscFunctionBegin;
1505d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1506d405a339SMatthew Knepley   if (!ts->setupcalled) {
1507d405a339SMatthew Knepley     ierr = TSSetUp(ts);                                                                                   CHKERRQ(ierr);
1508d405a339SMatthew Knepley   }
1509d405a339SMatthew Knepley 
1510d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);                                                        CHKERRQ(ierr);
1511d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);                                                                         CHKERRQ(ierr);
1512000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);                                                              CHKERRQ(ierr);
1513d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);                                                                        CHKERRQ(ierr);
1514d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);                                                          CHKERRQ(ierr);
1515d405a339SMatthew Knepley 
1516d405a339SMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);                                               CHKERRQ(ierr);
1517d405a339SMatthew Knepley   if ((opt == PETSC_TRUE) && !PetscPreLoadingOn) {
1518d405a339SMatthew Knepley     ierr = TSView(ts, PETSC_VIEWER_STDOUT_WORLD);                                                         CHKERRQ(ierr);
1519d405a339SMatthew Knepley   }
1520d763cef2SBarry Smith   PetscFunctionReturn(0);
1521d763cef2SBarry Smith }
1522d763cef2SBarry Smith 
15234a2ae208SSatish Balay #undef __FUNCT__
15244a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1525d763cef2SBarry Smith /*
1526d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1527d763cef2SBarry Smith */
152887828ca2SBarry Smith int TSMonitor(TS ts,int step,PetscReal ptime,Vec x)
1529d763cef2SBarry Smith {
1530d763cef2SBarry Smith   int i,ierr,n = ts->numbermonitors;
1531d763cef2SBarry Smith 
1532d763cef2SBarry Smith   PetscFunctionBegin;
1533d763cef2SBarry Smith   for (i=0; i<n; i++) {
1534142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1535d763cef2SBarry Smith   }
1536d763cef2SBarry Smith   PetscFunctionReturn(0);
1537d763cef2SBarry Smith }
1538d763cef2SBarry Smith 
1539d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1540d763cef2SBarry Smith 
15414a2ae208SSatish Balay #undef __FUNCT__
15424a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1543d763cef2SBarry Smith /*@C
1544d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1545d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1546d763cef2SBarry Smith 
1547d763cef2SBarry Smith    Collective on TS
1548d763cef2SBarry Smith 
1549d763cef2SBarry Smith    Input Parameters:
1550d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1551d763cef2SBarry Smith .  label - the title to put in the title bar
15527c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1553d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1554d763cef2SBarry Smith 
1555d763cef2SBarry Smith    Output Parameter:
1556d763cef2SBarry Smith .  draw - the drawing context
1557d763cef2SBarry Smith 
1558d763cef2SBarry Smith    Options Database Key:
1559d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1560d763cef2SBarry Smith 
1561d763cef2SBarry Smith    Notes:
1562b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1563d763cef2SBarry Smith 
1564d763cef2SBarry Smith    Level: intermediate
1565d763cef2SBarry Smith 
15667c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1567d763cef2SBarry Smith 
1568d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
15697c922b88SBarry Smith 
1570d763cef2SBarry Smith @*/
1571b0a32e0cSBarry Smith int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1572d763cef2SBarry Smith {
1573b0a32e0cSBarry Smith   PetscDraw win;
1574d763cef2SBarry Smith   int       ierr;
1575d763cef2SBarry Smith 
1576d763cef2SBarry Smith   PetscFunctionBegin;
1577b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1578b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1579b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1580b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1581d763cef2SBarry Smith 
1582b0a32e0cSBarry Smith   PetscLogObjectParent(*draw,win);
1583d763cef2SBarry Smith   PetscFunctionReturn(0);
1584d763cef2SBarry Smith }
1585d763cef2SBarry Smith 
15864a2ae208SSatish Balay #undef __FUNCT__
15874a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
158887828ca2SBarry Smith int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx)
1589d763cef2SBarry Smith {
1590b0a32e0cSBarry Smith   PetscDrawLG lg = (PetscDrawLG) monctx;
159187828ca2SBarry Smith   PetscReal      x,y = ptime;
1592d763cef2SBarry Smith   int         ierr;
1593d763cef2SBarry Smith 
1594d763cef2SBarry Smith   PetscFunctionBegin;
15957c922b88SBarry Smith   if (!monctx) {
15967c922b88SBarry Smith     MPI_Comm    comm;
1597b0a32e0cSBarry Smith     PetscViewer viewer;
15987c922b88SBarry Smith 
15997c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1600b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1601b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
16027c922b88SBarry Smith   }
16037c922b88SBarry Smith 
1604b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
160587828ca2SBarry Smith   x = (PetscReal)n;
1606b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1607d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1608b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1609d763cef2SBarry Smith   }
1610d763cef2SBarry Smith   PetscFunctionReturn(0);
1611d763cef2SBarry Smith }
1612d763cef2SBarry Smith 
16134a2ae208SSatish Balay #undef __FUNCT__
16144a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1615d763cef2SBarry Smith /*@C
1616d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1617d763cef2SBarry Smith    with TSLGMonitorCreate().
1618d763cef2SBarry Smith 
1619b0a32e0cSBarry Smith    Collective on PetscDrawLG
1620d763cef2SBarry Smith 
1621d763cef2SBarry Smith    Input Parameter:
1622d763cef2SBarry Smith .  draw - the drawing context
1623d763cef2SBarry Smith 
1624d763cef2SBarry Smith    Level: intermediate
1625d763cef2SBarry Smith 
1626d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1627d763cef2SBarry Smith 
1628d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1629d763cef2SBarry Smith @*/
1630b0a32e0cSBarry Smith int TSLGMonitorDestroy(PetscDrawLG drawlg)
1631d763cef2SBarry Smith {
1632b0a32e0cSBarry Smith   PetscDraw draw;
1633d763cef2SBarry Smith   int       ierr;
1634d763cef2SBarry Smith 
1635d763cef2SBarry Smith   PetscFunctionBegin;
1636b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1637b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1638b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1639d763cef2SBarry Smith   PetscFunctionReturn(0);
1640d763cef2SBarry Smith }
1641d763cef2SBarry Smith 
16424a2ae208SSatish Balay #undef __FUNCT__
16434a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1644d763cef2SBarry Smith /*@
1645d763cef2SBarry Smith    TSGetTime - Gets the current time.
1646d763cef2SBarry Smith 
1647d763cef2SBarry Smith    Not Collective
1648d763cef2SBarry Smith 
1649d763cef2SBarry Smith    Input Parameter:
1650d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1651d763cef2SBarry Smith 
1652d763cef2SBarry Smith    Output Parameter:
1653d763cef2SBarry Smith .  t  - the current time
1654d763cef2SBarry Smith 
1655d763cef2SBarry Smith    Contributed by: Matthew Knepley
1656d763cef2SBarry Smith 
1657d763cef2SBarry Smith    Level: beginner
1658d763cef2SBarry Smith 
1659d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1660d763cef2SBarry Smith 
1661d763cef2SBarry Smith .keywords: TS, get, time
1662d763cef2SBarry Smith @*/
166387828ca2SBarry Smith int TSGetTime(TS ts,PetscReal* t)
1664d763cef2SBarry Smith {
1665d763cef2SBarry Smith   PetscFunctionBegin;
1666d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1667d763cef2SBarry Smith   *t = ts->ptime;
1668d763cef2SBarry Smith   PetscFunctionReturn(0);
1669d763cef2SBarry Smith }
1670d763cef2SBarry Smith 
16714a2ae208SSatish Balay #undef __FUNCT__
16724a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1673d763cef2SBarry Smith /*@C
1674d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1675d763cef2SBarry Smith    TS options in the database.
1676d763cef2SBarry Smith 
1677d763cef2SBarry Smith    Collective on TS
1678d763cef2SBarry Smith 
1679d763cef2SBarry Smith    Input Parameter:
1680d763cef2SBarry Smith +  ts     - The TS context
1681d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1682d763cef2SBarry Smith 
1683d763cef2SBarry Smith    Notes:
1684d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1685d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1686d763cef2SBarry Smith    hyphen.
1687d763cef2SBarry Smith 
1688d763cef2SBarry Smith    Contributed by: Matthew Knepley
1689d763cef2SBarry Smith 
1690d763cef2SBarry Smith    Level: advanced
1691d763cef2SBarry Smith 
1692d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1693d763cef2SBarry Smith 
1694d763cef2SBarry Smith .seealso: TSSetFromOptions()
1695d763cef2SBarry Smith 
1696d763cef2SBarry Smith @*/
1697d763cef2SBarry Smith int TSSetOptionsPrefix(TS ts,char *prefix)
1698d763cef2SBarry Smith {
1699d763cef2SBarry Smith   int ierr;
1700d763cef2SBarry Smith 
1701d763cef2SBarry Smith   PetscFunctionBegin;
1702d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1703d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1704d763cef2SBarry Smith   switch(ts->problem_type) {
1705d763cef2SBarry Smith     case TS_NONLINEAR:
1706d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1707d763cef2SBarry Smith       break;
1708d763cef2SBarry Smith     case TS_LINEAR:
1709d763cef2SBarry Smith       ierr = SLESSetOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr);
1710d763cef2SBarry Smith       break;
1711d763cef2SBarry Smith   }
1712d763cef2SBarry Smith   PetscFunctionReturn(0);
1713d763cef2SBarry Smith }
1714d763cef2SBarry Smith 
1715d763cef2SBarry Smith 
17164a2ae208SSatish Balay #undef __FUNCT__
17174a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1718d763cef2SBarry Smith /*@C
1719d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1720d763cef2SBarry Smith    TS options in the database.
1721d763cef2SBarry Smith 
1722d763cef2SBarry Smith    Collective on TS
1723d763cef2SBarry Smith 
1724d763cef2SBarry Smith    Input Parameter:
1725d763cef2SBarry Smith +  ts     - The TS context
1726d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1727d763cef2SBarry Smith 
1728d763cef2SBarry Smith    Notes:
1729d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1730d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1731d763cef2SBarry Smith    hyphen.
1732d763cef2SBarry Smith 
1733d763cef2SBarry Smith    Contributed by: Matthew Knepley
1734d763cef2SBarry Smith 
1735d763cef2SBarry Smith    Level: advanced
1736d763cef2SBarry Smith 
1737d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1738d763cef2SBarry Smith 
1739d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1740d763cef2SBarry Smith 
1741d763cef2SBarry Smith @*/
1742d763cef2SBarry Smith int TSAppendOptionsPrefix(TS ts,char *prefix)
1743d763cef2SBarry Smith {
1744d763cef2SBarry Smith   int ierr;
1745d763cef2SBarry Smith 
1746d763cef2SBarry Smith   PetscFunctionBegin;
1747d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1748d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1749d763cef2SBarry Smith   switch(ts->problem_type) {
1750d763cef2SBarry Smith     case TS_NONLINEAR:
1751d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1752d763cef2SBarry Smith       break;
1753d763cef2SBarry Smith     case TS_LINEAR:
1754d763cef2SBarry Smith       ierr = SLESAppendOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr);
1755d763cef2SBarry Smith       break;
1756d763cef2SBarry Smith   }
1757d763cef2SBarry Smith   PetscFunctionReturn(0);
1758d763cef2SBarry Smith }
1759d763cef2SBarry Smith 
17604a2ae208SSatish Balay #undef __FUNCT__
17614a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1762d763cef2SBarry Smith /*@C
1763d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1764d763cef2SBarry Smith    TS options in the database.
1765d763cef2SBarry Smith 
1766d763cef2SBarry Smith    Not Collective
1767d763cef2SBarry Smith 
1768d763cef2SBarry Smith    Input Parameter:
1769d763cef2SBarry Smith .  ts - The TS context
1770d763cef2SBarry Smith 
1771d763cef2SBarry Smith    Output Parameter:
1772d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1773d763cef2SBarry Smith 
1774d763cef2SBarry Smith    Contributed by: Matthew Knepley
1775d763cef2SBarry Smith 
1776d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1777d763cef2SBarry Smith    sufficient length to hold the prefix.
1778d763cef2SBarry Smith 
1779d763cef2SBarry Smith    Level: intermediate
1780d763cef2SBarry Smith 
1781d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1782d763cef2SBarry Smith 
1783d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1784d763cef2SBarry Smith @*/
1785d763cef2SBarry Smith int TSGetOptionsPrefix(TS ts,char **prefix)
1786d763cef2SBarry Smith {
1787d763cef2SBarry Smith   int ierr;
1788d763cef2SBarry Smith 
1789d763cef2SBarry Smith   PetscFunctionBegin;
1790d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1791d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1792d763cef2SBarry Smith   PetscFunctionReturn(0);
1793d763cef2SBarry Smith }
1794d763cef2SBarry Smith 
17954a2ae208SSatish Balay #undef __FUNCT__
17964a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1797d763cef2SBarry Smith /*@C
1798d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1799d763cef2SBarry Smith 
1800d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1801d763cef2SBarry Smith 
1802d763cef2SBarry Smith    Input Parameter:
1803d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1804d763cef2SBarry Smith 
1805d763cef2SBarry Smith    Output Parameters:
1806d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1807d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1808d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1809d763cef2SBarry Smith 
1810d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1811d763cef2SBarry Smith 
1812d763cef2SBarry Smith    Contributed by: Matthew Knepley
1813d763cef2SBarry Smith 
1814d763cef2SBarry Smith    Level: intermediate
1815d763cef2SBarry Smith 
1816d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1817d763cef2SBarry Smith 
1818d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1819d763cef2SBarry Smith 
1820d763cef2SBarry Smith @*/
1821d763cef2SBarry Smith int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1822d763cef2SBarry Smith {
1823d763cef2SBarry Smith   PetscFunctionBegin;
1824d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1825d763cef2SBarry Smith   if (A)   *A = ts->A;
1826d763cef2SBarry Smith   if (M)   *M = ts->B;
1827d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1828d763cef2SBarry Smith   PetscFunctionReturn(0);
1829d763cef2SBarry Smith }
1830d763cef2SBarry Smith 
18314a2ae208SSatish Balay #undef __FUNCT__
18324a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1833d763cef2SBarry Smith /*@C
1834d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1835d763cef2SBarry Smith 
1836d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1837d763cef2SBarry Smith 
1838d763cef2SBarry Smith    Input Parameter:
1839d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1840d763cef2SBarry Smith 
1841d763cef2SBarry Smith    Output Parameters:
1842d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1843d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1844d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1845d763cef2SBarry Smith 
1846d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1847d763cef2SBarry Smith 
1848d763cef2SBarry Smith    Contributed by: Matthew Knepley
1849d763cef2SBarry Smith 
1850d763cef2SBarry Smith    Level: intermediate
1851d763cef2SBarry Smith 
1852d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1853d763cef2SBarry Smith 
1854d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1855d763cef2SBarry Smith @*/
1856d763cef2SBarry Smith int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1857d763cef2SBarry Smith {
1858d763cef2SBarry Smith   int ierr;
1859d763cef2SBarry Smith 
1860d763cef2SBarry Smith   PetscFunctionBegin;
1861d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1862d763cef2SBarry Smith   PetscFunctionReturn(0);
1863d763cef2SBarry Smith }
1864d763cef2SBarry Smith 
18651713a123SBarry Smith #undef __FUNCT__
18661713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
18671713a123SBarry Smith /*@C
18681713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
18691713a123SBarry Smith    VecView() for the solution at each timestep
18701713a123SBarry Smith 
18711713a123SBarry Smith    Collective on TS
18721713a123SBarry Smith 
18731713a123SBarry Smith    Input Parameters:
18741713a123SBarry Smith +  ts - the TS context
18751713a123SBarry Smith .  step - current time-step
1876142b95e3SSatish Balay .  ptime - current time
18771713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
18781713a123SBarry Smith 
18791713a123SBarry Smith    Level: intermediate
18801713a123SBarry Smith 
18811713a123SBarry Smith .keywords: TS,  vector, monitor, view
18821713a123SBarry Smith 
18831713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
18841713a123SBarry Smith @*/
1885142b95e3SSatish Balay int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
18861713a123SBarry Smith {
18871713a123SBarry Smith   int         ierr;
18881713a123SBarry Smith   PetscViewer viewer = (PetscViewer) dummy;
18891713a123SBarry Smith 
18901713a123SBarry Smith   PetscFunctionBegin;
18911713a123SBarry Smith   if (!viewer) {
18921713a123SBarry Smith     MPI_Comm comm;
18931713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
18941713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
18951713a123SBarry Smith   }
18961713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
18971713a123SBarry Smith   PetscFunctionReturn(0);
18981713a123SBarry Smith }
18991713a123SBarry Smith 
19001713a123SBarry Smith 
19011713a123SBarry Smith 
1902