xref: /petsc/src/ts/interface/ts.c (revision ef66eb6987ddfdf4e414d6b820cbc8d8d7d17bc2)
1*ef66eb69SBarry Smith /* $Id: ts.c,v 1.42 2001/08/06 21:18:08 bsmith Exp bsmith $ */
2e090d566SSatish Balay #include "src/ts/tsimpl.h"        /*I "petscts.h"  I*/
3d763cef2SBarry Smith 
44a2ae208SSatish Balay #undef __FUNCT__
54a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
6a7a1495cSBarry Smith /*@
78c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
8a7a1495cSBarry Smith       set with TSSetRHSJacobian().
9a7a1495cSBarry Smith 
10a7a1495cSBarry Smith    Collective on TS and Vec
11a7a1495cSBarry Smith 
12a7a1495cSBarry Smith    Input Parameters:
13a7a1495cSBarry Smith +  ts - the SNES context
14a7a1495cSBarry Smith .  t - current timestep
15a7a1495cSBarry Smith -  x - input vector
16a7a1495cSBarry Smith 
17a7a1495cSBarry Smith    Output Parameters:
18a7a1495cSBarry Smith +  A - Jacobian matrix
19a7a1495cSBarry Smith .  B - optional preconditioning matrix
20a7a1495cSBarry Smith -  flag - flag indicating matrix structure
21a7a1495cSBarry Smith 
22a7a1495cSBarry Smith    Notes:
23a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
24a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
25a7a1495cSBarry Smith 
26a7a1495cSBarry Smith    See SLESSetOperators() for important information about setting the
27a7a1495cSBarry Smith    flag parameter.
28a7a1495cSBarry Smith 
29a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
30a7a1495cSBarry Smith 
31a7a1495cSBarry Smith    Level: developer
32a7a1495cSBarry Smith 
33a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
34a7a1495cSBarry Smith 
35a7a1495cSBarry Smith .seealso:  TSSetRHSJacobian(), SLESSetOperators()
36a7a1495cSBarry Smith @*/
3787828ca2SBarry Smith int TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
38a7a1495cSBarry Smith {
39a7a1495cSBarry Smith   int ierr;
40a7a1495cSBarry Smith 
41a7a1495cSBarry Smith   PetscFunctionBegin;
42a7a1495cSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
43a7a1495cSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE);
44a7a1495cSBarry Smith   PetscCheckSameComm(ts,X);
45a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
4629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
47a7a1495cSBarry Smith   }
48*ef66eb69SBarry Smith   if (ts->rhsjacobian) {
49b0a32e0cSBarry Smith     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
50a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
51a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
52a7a1495cSBarry Smith     ierr = (*ts->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
53a7a1495cSBarry Smith     PetscStackPop;
54b0a32e0cSBarry Smith     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
55a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
56a7a1495cSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE);
57a7a1495cSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE);
58*ef66eb69SBarry Smith   } else {
59*ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60*ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61*ef66eb69SBarry Smith     if (*A != *B) {
62*ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63*ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64*ef66eb69SBarry Smith     }
65*ef66eb69SBarry Smith   }
66a7a1495cSBarry Smith   PetscFunctionReturn(0);
67a7a1495cSBarry Smith }
68a7a1495cSBarry Smith 
694a2ae208SSatish Balay #undef __FUNCT__
704a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
71d763cef2SBarry Smith /*
72d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
73d763cef2SBarry Smith 
74d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
75d763cef2SBarry Smith    this routine applies the matrix.
76d763cef2SBarry Smith */
7787828ca2SBarry Smith int TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
78d763cef2SBarry Smith {
79d763cef2SBarry Smith   int ierr;
80d763cef2SBarry Smith 
81d763cef2SBarry Smith   PetscFunctionBegin;
82d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
837c922b88SBarry Smith   PetscValidHeader(x);
847c922b88SBarry Smith   PetscValidHeader(y);
85d763cef2SBarry Smith 
86b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
87d763cef2SBarry Smith   if (ts->rhsfunction) {
88d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
89d763cef2SBarry Smith     ierr = (*ts->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
90d763cef2SBarry Smith     PetscStackPop;
917c922b88SBarry Smith   } else {
92d763cef2SBarry Smith     if (ts->rhsmatrix) { /* assemble matrix for this timestep */
93d763cef2SBarry Smith       MatStructure flg;
94d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
95d763cef2SBarry Smith       ierr = (*ts->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
96d763cef2SBarry Smith       PetscStackPop;
97d763cef2SBarry Smith     }
98d763cef2SBarry Smith     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
997c922b88SBarry Smith   }
100d763cef2SBarry Smith 
101d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
102d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
103b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
104d763cef2SBarry Smith 
105d763cef2SBarry Smith   PetscFunctionReturn(0);
106d763cef2SBarry Smith }
107d763cef2SBarry Smith 
1084a2ae208SSatish Balay #undef __FUNCT__
1094a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
110d763cef2SBarry Smith /*@C
111d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
112d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
113d763cef2SBarry Smith 
114d763cef2SBarry Smith     Collective on TS
115d763cef2SBarry Smith 
116d763cef2SBarry Smith     Input Parameters:
117d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
118d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
119d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
120d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
121d763cef2SBarry Smith 
122d763cef2SBarry Smith     Calling sequence of func:
12387828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
124d763cef2SBarry Smith 
125d763cef2SBarry Smith +   t - current timestep
126d763cef2SBarry Smith .   u - input vector
127d763cef2SBarry Smith .   F - function vector
128d763cef2SBarry Smith -   ctx - [optional] user-defined function context
129d763cef2SBarry Smith 
130d763cef2SBarry Smith     Important:
131d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
132d763cef2SBarry Smith 
133d763cef2SBarry Smith     Level: beginner
134d763cef2SBarry Smith 
135d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
136d763cef2SBarry Smith 
137d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
138d763cef2SBarry Smith @*/
13987828ca2SBarry Smith int TSSetRHSFunction(TS ts,int (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
140d763cef2SBarry Smith {
141d763cef2SBarry Smith   PetscFunctionBegin;
142d763cef2SBarry Smith 
143d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
144d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
14529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
146d763cef2SBarry Smith   }
147d763cef2SBarry Smith   ts->rhsfunction = f;
148d763cef2SBarry Smith   ts->funP        = ctx;
149d763cef2SBarry Smith   PetscFunctionReturn(0);
150d763cef2SBarry Smith }
151d763cef2SBarry Smith 
1524a2ae208SSatish Balay #undef __FUNCT__
1534a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
154d763cef2SBarry Smith /*@C
155d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
156d763cef2SBarry Smith    Also sets the location to store A.
157d763cef2SBarry Smith 
158d763cef2SBarry Smith    Collective on TS
159d763cef2SBarry Smith 
160d763cef2SBarry Smith    Input Parameters:
161d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
162d763cef2SBarry Smith .  A   - matrix
163d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
164d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
165d763cef2SBarry Smith          if A is not a function of t.
166d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
167d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
168d763cef2SBarry Smith 
169d763cef2SBarry Smith    Calling sequence of func:
17087828ca2SBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,int *flag,void *ctx);
171d763cef2SBarry Smith 
172d763cef2SBarry Smith +  t - current timestep
173d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
174d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
175d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
176d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
177d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
178d763cef2SBarry Smith 
179d763cef2SBarry Smith    Notes:
180d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
181d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
182d763cef2SBarry Smith 
183d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
184d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
185d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
186d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
187d763cef2SBarry Smith    throughout the global iterations.
188d763cef2SBarry Smith 
189d763cef2SBarry Smith    Important:
190d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
191d763cef2SBarry Smith 
192d763cef2SBarry Smith    Level: beginner
193d763cef2SBarry Smith 
194d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
195d763cef2SBarry Smith 
196d763cef2SBarry Smith .seealso: TSSetRHSFunction()
197d763cef2SBarry Smith @*/
19887828ca2SBarry Smith int TSSetRHSMatrix(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
199d763cef2SBarry Smith {
200d763cef2SBarry Smith   PetscFunctionBegin;
201d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
202184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
203184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
204184914b5SBarry Smith   PetscCheckSameComm(ts,A);
205184914b5SBarry Smith   PetscCheckSameComm(ts,B);
206d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
20729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
208d763cef2SBarry Smith   }
209d763cef2SBarry Smith 
210d763cef2SBarry Smith   ts->rhsmatrix = f;
211d763cef2SBarry Smith   ts->jacP      = ctx;
212d763cef2SBarry Smith   ts->A         = A;
213d763cef2SBarry Smith   ts->B         = B;
214d763cef2SBarry Smith 
215d763cef2SBarry Smith   PetscFunctionReturn(0);
216d763cef2SBarry Smith }
217d763cef2SBarry Smith 
2184a2ae208SSatish Balay #undef __FUNCT__
2194a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
220d763cef2SBarry Smith /*@C
221d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
222d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
223d763cef2SBarry Smith 
224d763cef2SBarry Smith    Collective on TS
225d763cef2SBarry Smith 
226d763cef2SBarry Smith    Input Parameters:
227d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
228d763cef2SBarry Smith .  A   - Jacobian matrix
229d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
230d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
231d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
232d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
233d763cef2SBarry Smith 
234d763cef2SBarry Smith    Calling sequence of func:
23587828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
236d763cef2SBarry Smith 
237d763cef2SBarry Smith +  t - current timestep
238d763cef2SBarry Smith .  u - input vector
239d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
240d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
241d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
242d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
243d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
244d763cef2SBarry Smith 
245d763cef2SBarry Smith    Notes:
246d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
247d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
248d763cef2SBarry Smith 
249d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
250d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
251d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
252d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
253d763cef2SBarry Smith    throughout the global iterations.
254d763cef2SBarry Smith 
255d763cef2SBarry Smith    Level: beginner
256d763cef2SBarry Smith 
257d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
258d763cef2SBarry Smith 
259d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
260d763cef2SBarry Smith           SNESDefaultComputeJacobianColor()
261d763cef2SBarry Smith 
262d763cef2SBarry Smith @*/
26387828ca2SBarry Smith int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
264d763cef2SBarry Smith {
265d763cef2SBarry Smith   PetscFunctionBegin;
266d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
267184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
268184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
269184914b5SBarry Smith   PetscCheckSameComm(ts,A);
270184914b5SBarry Smith   PetscCheckSameComm(ts,B);
271d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
27229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
273d763cef2SBarry Smith   }
274d763cef2SBarry Smith 
275d763cef2SBarry Smith   ts->rhsjacobian = f;
276d763cef2SBarry Smith   ts->jacP        = ctx;
277d763cef2SBarry Smith   ts->A           = A;
278d763cef2SBarry Smith   ts->B           = B;
279d763cef2SBarry Smith   PetscFunctionReturn(0);
280d763cef2SBarry Smith }
281d763cef2SBarry Smith 
2824a2ae208SSatish Balay #undef __FUNCT__
2834a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
284d763cef2SBarry Smith /*
285d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
286d763cef2SBarry Smith 
287d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
288d763cef2SBarry Smith    this routine applies the matrix.
289d763cef2SBarry Smith */
29087828ca2SBarry Smith int TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
291d763cef2SBarry Smith {
292d763cef2SBarry Smith   int ierr;
293d763cef2SBarry Smith 
294d763cef2SBarry Smith   PetscFunctionBegin;
295d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
296d763cef2SBarry Smith   PetscValidHeader(x);
297184914b5SBarry Smith   PetscCheckSameComm(ts,x);
298d763cef2SBarry Smith 
299d763cef2SBarry Smith   if (ts->rhsbc) {
300d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
301d763cef2SBarry Smith     ierr = (*ts->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
302d763cef2SBarry Smith     PetscStackPop;
303d763cef2SBarry Smith     PetscFunctionReturn(0);
304d763cef2SBarry Smith   }
305d763cef2SBarry Smith 
306d763cef2SBarry Smith   PetscFunctionReturn(0);
307d763cef2SBarry Smith }
308d763cef2SBarry Smith 
3094a2ae208SSatish Balay #undef __FUNCT__
3104a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
311d763cef2SBarry Smith /*@C
312d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
313d763cef2SBarry Smith     boundary conditions for the function F.
314d763cef2SBarry Smith 
315d763cef2SBarry Smith     Collective on TS
316d763cef2SBarry Smith 
317d763cef2SBarry Smith     Input Parameters:
318d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
319d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
320d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
321d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
322d763cef2SBarry Smith 
323d763cef2SBarry Smith     Calling sequence of func:
32487828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
325d763cef2SBarry Smith 
326d763cef2SBarry Smith +   t - current timestep
327d763cef2SBarry Smith .   F - function vector
328d763cef2SBarry Smith -   ctx - [optional] user-defined function context
329d763cef2SBarry Smith 
330d763cef2SBarry Smith     Level: intermediate
331d763cef2SBarry Smith 
332d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
333d763cef2SBarry Smith @*/
33487828ca2SBarry Smith int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,PetscReal,Vec,void*),void *ctx)
335d763cef2SBarry Smith {
336d763cef2SBarry Smith   PetscFunctionBegin;
337d763cef2SBarry Smith 
338d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
339d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
34029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
341d763cef2SBarry Smith   }
342d763cef2SBarry Smith   ts->rhsbc = f;
343d763cef2SBarry Smith   ts->bcP   = ctx;
344d763cef2SBarry Smith   PetscFunctionReturn(0);
345d763cef2SBarry Smith }
346d763cef2SBarry Smith 
3474a2ae208SSatish Balay #undef __FUNCT__
3484a2ae208SSatish Balay #define __FUNCT__ "TSView"
3497e2c5f70SBarry Smith /*@C
350d763cef2SBarry Smith     TSView - Prints the TS data structure.
351d763cef2SBarry Smith 
3524c49b128SBarry Smith     Collective on TS
353d763cef2SBarry Smith 
354d763cef2SBarry Smith     Input Parameters:
355d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
356d763cef2SBarry Smith -   viewer - visualization context
357d763cef2SBarry Smith 
358d763cef2SBarry Smith     Options Database Key:
359d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
360d763cef2SBarry Smith 
361d763cef2SBarry Smith     Notes:
362d763cef2SBarry Smith     The available visualization contexts include
363b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
364b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
365d763cef2SBarry Smith          output where only the first processor opens
366d763cef2SBarry Smith          the file.  All other processors send their
367d763cef2SBarry Smith          data to the first processor to print.
368d763cef2SBarry Smith 
369d763cef2SBarry Smith     The user can open an alternative visualization context with
370b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
371d763cef2SBarry Smith 
372d763cef2SBarry Smith     Level: beginner
373d763cef2SBarry Smith 
374d763cef2SBarry Smith .keywords: TS, timestep, view
375d763cef2SBarry Smith 
376b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
377d763cef2SBarry Smith @*/
378b0a32e0cSBarry Smith int TSView(TS ts,PetscViewer viewer)
379d763cef2SBarry Smith {
380d763cef2SBarry Smith   int        ierr;
381454a90a3SBarry Smith   char       *type;
3826831982aSBarry Smith   PetscTruth isascii,isstring;
383d763cef2SBarry Smith 
384d763cef2SBarry Smith   PetscFunctionBegin;
385d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
386b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
387b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
3886831982aSBarry Smith   PetscCheckSameComm(ts,viewer);
389fd16b177SBarry Smith 
390b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
391b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
3920f5bd95cSBarry Smith   if (isascii) {
393b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
394454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
395454a90a3SBarry Smith     if (type) {
396b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
397184914b5SBarry Smith     } else {
398b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
399184914b5SBarry Smith     }
400d763cef2SBarry Smith     if (ts->view) {
401b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
402d763cef2SBarry Smith       ierr = (*ts->view)(ts,viewer);CHKERRQ(ierr);
403b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
404d763cef2SBarry Smith     }
405b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
406b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
407d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
408b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
409d763cef2SBarry Smith     }
410b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
4110f5bd95cSBarry Smith   } else if (isstring) {
412454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
413b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
414d763cef2SBarry Smith   }
415b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
416d763cef2SBarry Smith   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
417d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
418b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
419d763cef2SBarry Smith   PetscFunctionReturn(0);
420d763cef2SBarry Smith }
421d763cef2SBarry Smith 
422d763cef2SBarry Smith 
4234a2ae208SSatish Balay #undef __FUNCT__
4244a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
425d763cef2SBarry Smith /*@C
426d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
427d763cef2SBarry Smith    the timesteppers.
428d763cef2SBarry Smith 
429d763cef2SBarry Smith    Collective on TS
430d763cef2SBarry Smith 
431d763cef2SBarry Smith    Input Parameters:
432d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
433d763cef2SBarry Smith -  usrP - optional user context
434d763cef2SBarry Smith 
435d763cef2SBarry Smith    Level: intermediate
436d763cef2SBarry Smith 
437d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
438d763cef2SBarry Smith 
439d763cef2SBarry Smith .seealso: TSGetApplicationContext()
440d763cef2SBarry Smith @*/
441d763cef2SBarry Smith int TSSetApplicationContext(TS ts,void *usrP)
442d763cef2SBarry Smith {
443d763cef2SBarry Smith   PetscFunctionBegin;
444d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
445d763cef2SBarry Smith   ts->user = usrP;
446d763cef2SBarry Smith   PetscFunctionReturn(0);
447d763cef2SBarry Smith }
448d763cef2SBarry Smith 
4494a2ae208SSatish Balay #undef __FUNCT__
4504a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
451d763cef2SBarry Smith /*@C
452d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
453d763cef2SBarry Smith     timestepper.
454d763cef2SBarry Smith 
455d763cef2SBarry Smith     Not Collective
456d763cef2SBarry Smith 
457d763cef2SBarry Smith     Input Parameter:
458d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
459d763cef2SBarry Smith 
460d763cef2SBarry Smith     Output Parameter:
461d763cef2SBarry Smith .   usrP - user context
462d763cef2SBarry Smith 
463d763cef2SBarry Smith     Level: intermediate
464d763cef2SBarry Smith 
465d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
466d763cef2SBarry Smith 
467d763cef2SBarry Smith .seealso: TSSetApplicationContext()
468d763cef2SBarry Smith @*/
469d763cef2SBarry Smith int TSGetApplicationContext(TS ts,void **usrP)
470d763cef2SBarry Smith {
471d763cef2SBarry Smith   PetscFunctionBegin;
472d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
473d763cef2SBarry Smith   *usrP = ts->user;
474d763cef2SBarry Smith   PetscFunctionReturn(0);
475d763cef2SBarry Smith }
476d763cef2SBarry Smith 
4774a2ae208SSatish Balay #undef __FUNCT__
4784a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
479d763cef2SBarry Smith /*@
480d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
481d763cef2SBarry Smith 
482d763cef2SBarry Smith    Not Collective
483d763cef2SBarry Smith 
484d763cef2SBarry Smith    Input Parameter:
485d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
486d763cef2SBarry Smith 
487d763cef2SBarry Smith    Output Parameter:
488d763cef2SBarry Smith .  iter - number steps so far
489d763cef2SBarry Smith 
490d763cef2SBarry Smith    Level: intermediate
491d763cef2SBarry Smith 
492d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
493d763cef2SBarry Smith @*/
494d763cef2SBarry Smith int TSGetTimeStepNumber(TS ts,int* iter)
495d763cef2SBarry Smith {
496d763cef2SBarry Smith   PetscFunctionBegin;
497d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
498d763cef2SBarry Smith   *iter = ts->steps;
499d763cef2SBarry Smith   PetscFunctionReturn(0);
500d763cef2SBarry Smith }
501d763cef2SBarry Smith 
5024a2ae208SSatish Balay #undef __FUNCT__
5034a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
504d763cef2SBarry Smith /*@
505d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
506d763cef2SBarry Smith    as well as the initial time.
507d763cef2SBarry Smith 
508d763cef2SBarry Smith    Collective on TS
509d763cef2SBarry Smith 
510d763cef2SBarry Smith    Input Parameters:
511d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
512d763cef2SBarry Smith .  initial_time - the initial time
513d763cef2SBarry Smith -  time_step - the size of the timestep
514d763cef2SBarry Smith 
515d763cef2SBarry Smith    Level: intermediate
516d763cef2SBarry Smith 
517d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
518d763cef2SBarry Smith 
519d763cef2SBarry Smith .keywords: TS, set, initial, timestep
520d763cef2SBarry Smith @*/
52187828ca2SBarry Smith int TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
522d763cef2SBarry Smith {
523d763cef2SBarry Smith   PetscFunctionBegin;
524d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
525d763cef2SBarry Smith   ts->time_step         = time_step;
526d763cef2SBarry Smith   ts->initial_time_step = time_step;
527d763cef2SBarry Smith   ts->ptime             = initial_time;
528d763cef2SBarry Smith   PetscFunctionReturn(0);
529d763cef2SBarry Smith }
530d763cef2SBarry Smith 
5314a2ae208SSatish Balay #undef __FUNCT__
5324a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
533d763cef2SBarry Smith /*@
534d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
535d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
536d763cef2SBarry Smith 
537d763cef2SBarry Smith    Collective on TS
538d763cef2SBarry Smith 
539d763cef2SBarry Smith    Input Parameters:
540d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
541d763cef2SBarry Smith -  time_step - the size of the timestep
542d763cef2SBarry Smith 
543d763cef2SBarry Smith    Level: intermediate
544d763cef2SBarry Smith 
545d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
546d763cef2SBarry Smith 
547d763cef2SBarry Smith .keywords: TS, set, timestep
548d763cef2SBarry Smith @*/
54987828ca2SBarry Smith int TSSetTimeStep(TS ts,PetscReal time_step)
550d763cef2SBarry Smith {
551d763cef2SBarry Smith   PetscFunctionBegin;
552d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
553d763cef2SBarry Smith   ts->time_step = time_step;
554d763cef2SBarry Smith   PetscFunctionReturn(0);
555d763cef2SBarry Smith }
556d763cef2SBarry Smith 
5574a2ae208SSatish Balay #undef __FUNCT__
5584a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
559d763cef2SBarry Smith /*@
560d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
561d763cef2SBarry Smith 
562d763cef2SBarry Smith    Not Collective
563d763cef2SBarry Smith 
564d763cef2SBarry Smith    Input Parameter:
565d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
566d763cef2SBarry Smith 
567d763cef2SBarry Smith    Output Parameter:
568d763cef2SBarry Smith .  dt - the current timestep size
569d763cef2SBarry Smith 
570d763cef2SBarry Smith    Level: intermediate
571d763cef2SBarry Smith 
572d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
573d763cef2SBarry Smith 
574d763cef2SBarry Smith .keywords: TS, get, timestep
575d763cef2SBarry Smith @*/
57687828ca2SBarry Smith int TSGetTimeStep(TS ts,PetscReal* dt)
577d763cef2SBarry Smith {
578d763cef2SBarry Smith   PetscFunctionBegin;
579d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
580d763cef2SBarry Smith   *dt = ts->time_step;
581d763cef2SBarry Smith   PetscFunctionReturn(0);
582d763cef2SBarry Smith }
583d763cef2SBarry Smith 
5844a2ae208SSatish Balay #undef __FUNCT__
5854a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
586d763cef2SBarry Smith /*@C
587d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
588d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
589d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
590d763cef2SBarry Smith    the solution at the next timestep has been calculated.
591d763cef2SBarry Smith 
592d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
593d763cef2SBarry Smith 
594d763cef2SBarry Smith    Input Parameter:
595d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
596d763cef2SBarry Smith 
597d763cef2SBarry Smith    Output Parameter:
598d763cef2SBarry Smith .  v - the vector containing the solution
599d763cef2SBarry Smith 
600d763cef2SBarry Smith    Level: intermediate
601d763cef2SBarry Smith 
602d763cef2SBarry Smith .seealso: TSGetTimeStep()
603d763cef2SBarry Smith 
604d763cef2SBarry Smith .keywords: TS, timestep, get, solution
605d763cef2SBarry Smith @*/
606d763cef2SBarry Smith int TSGetSolution(TS ts,Vec *v)
607d763cef2SBarry Smith {
608d763cef2SBarry Smith   PetscFunctionBegin;
609d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
610d763cef2SBarry Smith   *v = ts->vec_sol_always;
611d763cef2SBarry Smith   PetscFunctionReturn(0);
612d763cef2SBarry Smith }
613d763cef2SBarry Smith 
6144a2ae208SSatish Balay #undef __FUNCT__
6154a2ae208SSatish Balay #define __FUNCT__ "TSPublish_Petsc"
616454a90a3SBarry Smith static int TSPublish_Petsc(PetscObject obj)
617d763cef2SBarry Smith {
618aa482453SBarry Smith #if defined(PETSC_HAVE_AMS)
619454a90a3SBarry Smith   TS   v = (TS) obj;
620d763cef2SBarry Smith   int  ierr;
62143d6d2cbSBarry Smith #endif
622d763cef2SBarry Smith 
623d763cef2SBarry Smith   PetscFunctionBegin;
624d763cef2SBarry Smith 
62543d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS)
626d763cef2SBarry Smith   /* if it is already published then return */
627d763cef2SBarry Smith   if (v->amem >=0) PetscFunctionReturn(0);
628d763cef2SBarry Smith 
629454a90a3SBarry Smith   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
630d763cef2SBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ,
631d763cef2SBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
632d763cef2SBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ,
633d763cef2SBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
634d763cef2SBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1,
635d763cef2SBarry Smith                                AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
636454a90a3SBarry Smith   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
637d763cef2SBarry Smith #endif
638d763cef2SBarry Smith   PetscFunctionReturn(0);
639d763cef2SBarry Smith }
640d763cef2SBarry Smith 
641d763cef2SBarry Smith /* -----------------------------------------------------------*/
642d763cef2SBarry Smith 
6434a2ae208SSatish Balay #undef __FUNCT__
6444a2ae208SSatish Balay #define __FUNCT__ "TSCreate"
645d763cef2SBarry Smith /*@C
646d763cef2SBarry Smith    TSCreate - Creates a timestepper context.
647d763cef2SBarry Smith 
648d763cef2SBarry Smith    Collective on MPI_Comm
649d763cef2SBarry Smith 
650d763cef2SBarry Smith    Input Parameter:
651d763cef2SBarry Smith +  comm - MPI communicator
652d763cef2SBarry Smith -  type - One of  TS_LINEAR,TS_NONLINEAR
653d763cef2SBarry Smith    where these types refer to problems of the forms
654d763cef2SBarry Smith .vb
655d763cef2SBarry Smith          U_t = A U
656d763cef2SBarry Smith          U_t = A(t) U
657d763cef2SBarry Smith          U_t = F(t,U)
658d763cef2SBarry Smith .ve
659d763cef2SBarry Smith 
660d763cef2SBarry Smith    Output Parameter:
661d763cef2SBarry Smith .  outts - the new TS context
662d763cef2SBarry Smith 
663d763cef2SBarry Smith    Level: beginner
664d763cef2SBarry Smith 
665d763cef2SBarry Smith .keywords: TS, timestep, create, context
666d763cef2SBarry Smith 
667435da068SBarry Smith .seealso: TSSetUp(), TSStep(), TSDestroy(), TSProblemType, TS
668d763cef2SBarry Smith @*/
669d763cef2SBarry Smith int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts)
670d763cef2SBarry Smith {
671d763cef2SBarry Smith   TS ts;
672d763cef2SBarry Smith 
673d763cef2SBarry Smith   PetscFunctionBegin;
674d763cef2SBarry Smith   *outts = 0;
675d763cef2SBarry Smith   PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView);
676b0a32e0cSBarry Smith   PetscLogObjectCreate(ts);
677d763cef2SBarry Smith   ts->bops->publish     = TSPublish_Petsc;
678d763cef2SBarry Smith   ts->max_steps         = 5000;
679d763cef2SBarry Smith   ts->max_time          = 5.0;
680d763cef2SBarry Smith   ts->time_step         = .1;
681d763cef2SBarry Smith   ts->initial_time_step = ts->time_step;
682d763cef2SBarry Smith   ts->steps             = 0;
683d763cef2SBarry Smith   ts->ptime             = 0.0;
684d763cef2SBarry Smith   ts->data              = 0;
685d763cef2SBarry Smith   ts->view              = 0;
686d763cef2SBarry Smith   ts->setupcalled       = 0;
687d763cef2SBarry Smith   ts->problem_type      = problemtype;
688d763cef2SBarry Smith   ts->numbermonitors    = 0;
689d763cef2SBarry Smith   ts->linear_its        = 0;
690d763cef2SBarry Smith   ts->nonlinear_its     = 0;
691d763cef2SBarry Smith 
692d763cef2SBarry Smith   *outts = ts;
693d763cef2SBarry Smith   PetscFunctionReturn(0);
694d763cef2SBarry Smith }
695d763cef2SBarry Smith 
696d763cef2SBarry Smith /* ----- Routines to initialize and destroy a timestepper ---- */
697d763cef2SBarry Smith 
6984a2ae208SSatish Balay #undef __FUNCT__
6994a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
700d763cef2SBarry Smith /*@
701d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
702d763cef2SBarry Smith    of a timestepper.
703d763cef2SBarry Smith 
704d763cef2SBarry Smith    Collective on TS
705d763cef2SBarry Smith 
706d763cef2SBarry Smith    Input Parameter:
707d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
708d763cef2SBarry Smith 
709d763cef2SBarry Smith    Notes:
710d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
711d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
712d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
713d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
714d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
715d763cef2SBarry Smith 
716d763cef2SBarry Smith    Level: advanced
717d763cef2SBarry Smith 
718d763cef2SBarry Smith .keywords: TS, timestep, setup
719d763cef2SBarry Smith 
720d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
721d763cef2SBarry Smith @*/
722d763cef2SBarry Smith int TSSetUp(TS ts)
723d763cef2SBarry Smith {
724d763cef2SBarry Smith   int ierr;
725d763cef2SBarry Smith 
726d763cef2SBarry Smith   PetscFunctionBegin;
727d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
72829bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
729d763cef2SBarry Smith   if (!ts->type_name) {
730d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
731d763cef2SBarry Smith   }
732d763cef2SBarry Smith   ierr = (*ts->setup)(ts);CHKERRQ(ierr);
733d763cef2SBarry Smith   ts->setupcalled = 1;
734d763cef2SBarry Smith   PetscFunctionReturn(0);
735d763cef2SBarry Smith }
736d763cef2SBarry Smith 
7374a2ae208SSatish Balay #undef __FUNCT__
7384a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
739d763cef2SBarry Smith /*@C
740d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
741d763cef2SBarry Smith    with TSCreate().
742d763cef2SBarry Smith 
743d763cef2SBarry Smith    Collective on TS
744d763cef2SBarry Smith 
745d763cef2SBarry Smith    Input Parameter:
746d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
747d763cef2SBarry Smith 
748d763cef2SBarry Smith    Level: beginner
749d763cef2SBarry Smith 
750d763cef2SBarry Smith .keywords: TS, timestepper, destroy
751d763cef2SBarry Smith 
752d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
753d763cef2SBarry Smith @*/
754d763cef2SBarry Smith int TSDestroy(TS ts)
755d763cef2SBarry Smith {
756329f5518SBarry Smith   int ierr,i;
757d763cef2SBarry Smith 
758d763cef2SBarry Smith   PetscFunctionBegin;
759d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
760d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
761d763cef2SBarry Smith 
762be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
7630f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
764be0abb6dSBarry Smith 
765d763cef2SBarry Smith   if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);}
766d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
767d763cef2SBarry Smith   ierr = (*(ts)->destroy)(ts);CHKERRQ(ierr);
768329f5518SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
769329f5518SBarry Smith     if (ts->mdestroy[i]) {
770329f5518SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
771329f5518SBarry Smith     }
772329f5518SBarry Smith   }
773b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)ts);
774d763cef2SBarry Smith   PetscHeaderDestroy((PetscObject)ts);
775d763cef2SBarry Smith   PetscFunctionReturn(0);
776d763cef2SBarry Smith }
777d763cef2SBarry Smith 
7784a2ae208SSatish Balay #undef __FUNCT__
7794a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
780d763cef2SBarry Smith /*@C
781d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
782d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
783d763cef2SBarry Smith 
784d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
785d763cef2SBarry Smith 
786d763cef2SBarry Smith    Input Parameter:
787d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
788d763cef2SBarry Smith 
789d763cef2SBarry Smith    Output Parameter:
790d763cef2SBarry Smith .  snes - the nonlinear solver context
791d763cef2SBarry Smith 
792d763cef2SBarry Smith    Notes:
793d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
794d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
795d763cef2SBarry Smith    SLES, KSP, and PC contexts as well.
796d763cef2SBarry Smith 
797d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
798d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
799d763cef2SBarry Smith 
800d763cef2SBarry Smith    Level: beginner
801d763cef2SBarry Smith 
802d763cef2SBarry Smith .keywords: timestep, get, SNES
803d763cef2SBarry Smith @*/
804d763cef2SBarry Smith int TSGetSNES(TS ts,SNES *snes)
805d763cef2SBarry Smith {
806d763cef2SBarry Smith   PetscFunctionBegin;
807d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
80829bbc08cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetSLES()");
809d763cef2SBarry Smith   *snes = ts->snes;
810d763cef2SBarry Smith   PetscFunctionReturn(0);
811d763cef2SBarry Smith }
812d763cef2SBarry Smith 
8134a2ae208SSatish Balay #undef __FUNCT__
8144a2ae208SSatish Balay #define __FUNCT__ "TSGetSLES"
815d763cef2SBarry Smith /*@C
816d763cef2SBarry Smith    TSGetSLES - Returns the SLES (linear solver) associated with
817d763cef2SBarry Smith    a TS (timestepper) context.
818d763cef2SBarry Smith 
819d763cef2SBarry Smith    Not Collective, but SLES is parallel if TS is parallel
820d763cef2SBarry Smith 
821d763cef2SBarry Smith    Input Parameter:
822d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
823d763cef2SBarry Smith 
824d763cef2SBarry Smith    Output Parameter:
825d763cef2SBarry Smith .  sles - the nonlinear solver context
826d763cef2SBarry Smith 
827d763cef2SBarry Smith    Notes:
828d763cef2SBarry Smith    The user can then directly manipulate the SLES context to set various
829d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
830d763cef2SBarry Smith    KSP and PC contexts as well.
831d763cef2SBarry Smith 
832d763cef2SBarry Smith    TSGetSLES() does not work for integrators that do not use SLES;
833d763cef2SBarry Smith    in this case TSGetSLES() returns PETSC_NULL in sles.
834d763cef2SBarry Smith 
835d763cef2SBarry Smith    Level: beginner
836d763cef2SBarry Smith 
837d763cef2SBarry Smith .keywords: timestep, get, SLES
838d763cef2SBarry Smith @*/
839d763cef2SBarry Smith int TSGetSLES(TS ts,SLES *sles)
840d763cef2SBarry Smith {
841d763cef2SBarry Smith   PetscFunctionBegin;
842d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
84329bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
844d763cef2SBarry Smith   *sles = ts->sles;
845d763cef2SBarry Smith   PetscFunctionReturn(0);
846d763cef2SBarry Smith }
847d763cef2SBarry Smith 
848d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
849d763cef2SBarry Smith 
8504a2ae208SSatish Balay #undef __FUNCT__
8514a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
852d763cef2SBarry Smith /*@
853d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
854d763cef2SBarry Smith    maximum time for iteration.
855d763cef2SBarry Smith 
856d763cef2SBarry Smith    Collective on TS
857d763cef2SBarry Smith 
858d763cef2SBarry Smith    Input Parameters:
859d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
860d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
861d763cef2SBarry Smith -  maxtime - final time to iterate to
862d763cef2SBarry Smith 
863d763cef2SBarry Smith    Options Database Keys:
864d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
865d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
866d763cef2SBarry Smith 
867d763cef2SBarry Smith    Notes:
868d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
869d763cef2SBarry Smith 
870d763cef2SBarry Smith    Level: intermediate
871d763cef2SBarry Smith 
872d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
873d763cef2SBarry Smith @*/
87487828ca2SBarry Smith int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime)
875d763cef2SBarry Smith {
876d763cef2SBarry Smith   PetscFunctionBegin;
877d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
878d763cef2SBarry Smith   ts->max_steps = maxsteps;
879d763cef2SBarry Smith   ts->max_time  = maxtime;
880d763cef2SBarry Smith   PetscFunctionReturn(0);
881d763cef2SBarry Smith }
882d763cef2SBarry Smith 
8834a2ae208SSatish Balay #undef __FUNCT__
8844a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
885d763cef2SBarry Smith /*@
886d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
887d763cef2SBarry Smith    for use by the TS routines.
888d763cef2SBarry Smith 
889d763cef2SBarry Smith    Collective on TS and Vec
890d763cef2SBarry Smith 
891d763cef2SBarry Smith    Input Parameters:
892d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
893d763cef2SBarry Smith -  x - the solution vector
894d763cef2SBarry Smith 
895d763cef2SBarry Smith    Level: beginner
896d763cef2SBarry Smith 
897d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
898d763cef2SBarry Smith @*/
899d763cef2SBarry Smith int TSSetSolution(TS ts,Vec x)
900d763cef2SBarry Smith {
901d763cef2SBarry Smith   PetscFunctionBegin;
902d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
903d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
904d763cef2SBarry Smith   PetscFunctionReturn(0);
905d763cef2SBarry Smith }
906d763cef2SBarry Smith 
907d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
908d763cef2SBarry Smith 
9094a2ae208SSatish Balay #undef __FUNCT__
9104a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
911d763cef2SBarry Smith /*@C
912d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
913d763cef2SBarry Smith    timestep to display the iteration's  progress.
914d763cef2SBarry Smith 
915d763cef2SBarry Smith    Collective on TS
916d763cef2SBarry Smith 
917d763cef2SBarry Smith    Input Parameters:
918d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
919d763cef2SBarry Smith .  func - monitoring routine
920329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
921b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
922b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
923b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
924d763cef2SBarry Smith 
925d763cef2SBarry Smith    Calling sequence of func:
92687828ca2SBarry Smith $    int func(TS ts,int steps,PetscReal time,Vec x,void *mctx)
927d763cef2SBarry Smith 
928d763cef2SBarry Smith +    ts - the TS context
929d763cef2SBarry Smith .    steps - iteration number
9301f06c33eSBarry Smith .    time - current time
931d763cef2SBarry Smith .    x - current iterate
932d763cef2SBarry Smith -    mctx - [optional] monitoring context
933d763cef2SBarry Smith 
934d763cef2SBarry Smith    Notes:
935d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
936d763cef2SBarry Smith    already has been loaded.
937d763cef2SBarry Smith 
938d763cef2SBarry Smith    Level: intermediate
939d763cef2SBarry Smith 
940d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
941d763cef2SBarry Smith 
942d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
943d763cef2SBarry Smith @*/
94487828ca2SBarry Smith int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*))
945d763cef2SBarry Smith {
946d763cef2SBarry Smith   PetscFunctionBegin;
947d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
948d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
94929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
950d763cef2SBarry Smith   }
951d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
952329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
953d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
954d763cef2SBarry Smith   PetscFunctionReturn(0);
955d763cef2SBarry Smith }
956d763cef2SBarry Smith 
9574a2ae208SSatish Balay #undef __FUNCT__
9584a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
959d763cef2SBarry Smith /*@C
960d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
961d763cef2SBarry Smith 
962d763cef2SBarry Smith    Collective on TS
963d763cef2SBarry Smith 
964d763cef2SBarry Smith    Input Parameters:
965d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
966d763cef2SBarry Smith 
967d763cef2SBarry Smith    Notes:
968d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
969d763cef2SBarry Smith 
970d763cef2SBarry Smith    Level: intermediate
971d763cef2SBarry Smith 
972d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
973d763cef2SBarry Smith 
974d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
975d763cef2SBarry Smith @*/
976d763cef2SBarry Smith int TSClearMonitor(TS ts)
977d763cef2SBarry Smith {
978d763cef2SBarry Smith   PetscFunctionBegin;
979d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
980d763cef2SBarry Smith   ts->numbermonitors = 0;
981d763cef2SBarry Smith   PetscFunctionReturn(0);
982d763cef2SBarry Smith }
983d763cef2SBarry Smith 
9844a2ae208SSatish Balay #undef __FUNCT__
9854a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
98687828ca2SBarry Smith int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
987d763cef2SBarry Smith {
988d132466eSBarry Smith   int ierr;
989d132466eSBarry Smith 
990d763cef2SBarry Smith   PetscFunctionBegin;
991142b95e3SSatish Balay   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
992d763cef2SBarry Smith   PetscFunctionReturn(0);
993d763cef2SBarry Smith }
994d763cef2SBarry Smith 
9954a2ae208SSatish Balay #undef __FUNCT__
9964a2ae208SSatish Balay #define __FUNCT__ "TSStep"
997d763cef2SBarry Smith /*@
998d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
999d763cef2SBarry Smith 
1000d763cef2SBarry Smith    Collective on TS
1001d763cef2SBarry Smith 
1002d763cef2SBarry Smith    Input Parameter:
1003d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1004d763cef2SBarry Smith 
1005d763cef2SBarry Smith    Output Parameters:
1006d763cef2SBarry Smith +  steps - number of iterations until termination
1007142b95e3SSatish Balay -  ptime - time until termination
1008d763cef2SBarry Smith 
1009d763cef2SBarry Smith    Level: beginner
1010d763cef2SBarry Smith 
1011d763cef2SBarry Smith .keywords: TS, timestep, solve
1012d763cef2SBarry Smith 
1013d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1014d763cef2SBarry Smith @*/
101587828ca2SBarry Smith int TSStep(TS ts,int *steps,PetscReal *ptime)
1016d763cef2SBarry Smith {
1017f1af5d2fSBarry Smith   int        ierr;
1018f1af5d2fSBarry Smith   PetscTruth flg;
1019d763cef2SBarry Smith 
1020d763cef2SBarry Smith   PetscFunctionBegin;
1021d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1022d763cef2SBarry Smith   if (!ts->setupcalled) {ierr = TSSetUp(ts);CHKERRQ(ierr);}
1023b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1024142b95e3SSatish Balay   ierr = (*ts->step)(ts,steps,ptime);CHKERRQ(ierr);
1025b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1026b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(ts->prefix,"-ts_view",&flg);CHKERRQ(ierr);
1027b5758dffSBarry Smith   if (flg  && !PetscPreLoadingOn) {ierr = TSView(ts,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
1028d763cef2SBarry Smith   PetscFunctionReturn(0);
1029d763cef2SBarry Smith }
1030d763cef2SBarry Smith 
10314a2ae208SSatish Balay #undef __FUNCT__
10324a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1033d763cef2SBarry Smith /*
1034d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1035d763cef2SBarry Smith */
103687828ca2SBarry Smith int TSMonitor(TS ts,int step,PetscReal ptime,Vec x)
1037d763cef2SBarry Smith {
1038d763cef2SBarry Smith   int i,ierr,n = ts->numbermonitors;
1039d763cef2SBarry Smith 
1040d763cef2SBarry Smith   PetscFunctionBegin;
1041d763cef2SBarry Smith   for (i=0; i<n; i++) {
1042142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1043d763cef2SBarry Smith   }
1044d763cef2SBarry Smith   PetscFunctionReturn(0);
1045d763cef2SBarry Smith }
1046d763cef2SBarry Smith 
1047d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1048d763cef2SBarry Smith 
10494a2ae208SSatish Balay #undef __FUNCT__
10504a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1051d763cef2SBarry Smith /*@C
1052d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1053d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1054d763cef2SBarry Smith 
1055d763cef2SBarry Smith    Collective on TS
1056d763cef2SBarry Smith 
1057d763cef2SBarry Smith    Input Parameters:
1058d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1059d763cef2SBarry Smith .  label - the title to put in the title bar
10607c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1061d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1062d763cef2SBarry Smith 
1063d763cef2SBarry Smith    Output Parameter:
1064d763cef2SBarry Smith .  draw - the drawing context
1065d763cef2SBarry Smith 
1066d763cef2SBarry Smith    Options Database Key:
1067d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1068d763cef2SBarry Smith 
1069d763cef2SBarry Smith    Notes:
1070b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1071d763cef2SBarry Smith 
1072d763cef2SBarry Smith    Level: intermediate
1073d763cef2SBarry Smith 
10747c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1075d763cef2SBarry Smith 
1076d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
10777c922b88SBarry Smith 
1078d763cef2SBarry Smith @*/
1079b0a32e0cSBarry Smith int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1080d763cef2SBarry Smith {
1081b0a32e0cSBarry Smith   PetscDraw win;
1082d763cef2SBarry Smith   int       ierr;
1083d763cef2SBarry Smith 
1084d763cef2SBarry Smith   PetscFunctionBegin;
1085b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1086b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1087b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1088b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1089d763cef2SBarry Smith 
1090b0a32e0cSBarry Smith   PetscLogObjectParent(*draw,win);
1091d763cef2SBarry Smith   PetscFunctionReturn(0);
1092d763cef2SBarry Smith }
1093d763cef2SBarry Smith 
10944a2ae208SSatish Balay #undef __FUNCT__
10954a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
109687828ca2SBarry Smith int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx)
1097d763cef2SBarry Smith {
1098b0a32e0cSBarry Smith   PetscDrawLG lg = (PetscDrawLG) monctx;
109987828ca2SBarry Smith   PetscReal      x,y = ptime;
1100d763cef2SBarry Smith   int         ierr;
1101d763cef2SBarry Smith 
1102d763cef2SBarry Smith   PetscFunctionBegin;
11037c922b88SBarry Smith   if (!monctx) {
11047c922b88SBarry Smith     MPI_Comm    comm;
1105b0a32e0cSBarry Smith     PetscViewer viewer;
11067c922b88SBarry Smith 
11077c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1108b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1109b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
11107c922b88SBarry Smith   }
11117c922b88SBarry Smith 
1112b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
111387828ca2SBarry Smith   x = (PetscReal)n;
1114b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1115d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1116b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1117d763cef2SBarry Smith   }
1118d763cef2SBarry Smith   PetscFunctionReturn(0);
1119d763cef2SBarry Smith }
1120d763cef2SBarry Smith 
11214a2ae208SSatish Balay #undef __FUNCT__
11224a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1123d763cef2SBarry Smith /*@C
1124d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1125d763cef2SBarry Smith    with TSLGMonitorCreate().
1126d763cef2SBarry Smith 
1127b0a32e0cSBarry Smith    Collective on PetscDrawLG
1128d763cef2SBarry Smith 
1129d763cef2SBarry Smith    Input Parameter:
1130d763cef2SBarry Smith .  draw - the drawing context
1131d763cef2SBarry Smith 
1132d763cef2SBarry Smith    Level: intermediate
1133d763cef2SBarry Smith 
1134d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1135d763cef2SBarry Smith 
1136d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1137d763cef2SBarry Smith @*/
1138b0a32e0cSBarry Smith int TSLGMonitorDestroy(PetscDrawLG drawlg)
1139d763cef2SBarry Smith {
1140b0a32e0cSBarry Smith   PetscDraw draw;
1141d763cef2SBarry Smith   int       ierr;
1142d763cef2SBarry Smith 
1143d763cef2SBarry Smith   PetscFunctionBegin;
1144b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1145b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1146b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1147d763cef2SBarry Smith   PetscFunctionReturn(0);
1148d763cef2SBarry Smith }
1149d763cef2SBarry Smith 
11504a2ae208SSatish Balay #undef __FUNCT__
11514a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1152d763cef2SBarry Smith /*@
1153d763cef2SBarry Smith    TSGetTime - Gets the current time.
1154d763cef2SBarry Smith 
1155d763cef2SBarry Smith    Not Collective
1156d763cef2SBarry Smith 
1157d763cef2SBarry Smith    Input Parameter:
1158d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1159d763cef2SBarry Smith 
1160d763cef2SBarry Smith    Output Parameter:
1161d763cef2SBarry Smith .  t  - the current time
1162d763cef2SBarry Smith 
1163d763cef2SBarry Smith    Contributed by: Matthew Knepley
1164d763cef2SBarry Smith 
1165d763cef2SBarry Smith    Level: beginner
1166d763cef2SBarry Smith 
1167d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1168d763cef2SBarry Smith 
1169d763cef2SBarry Smith .keywords: TS, get, time
1170d763cef2SBarry Smith @*/
117187828ca2SBarry Smith int TSGetTime(TS ts,PetscReal* t)
1172d763cef2SBarry Smith {
1173d763cef2SBarry Smith   PetscFunctionBegin;
1174d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1175d763cef2SBarry Smith   *t = ts->ptime;
1176d763cef2SBarry Smith   PetscFunctionReturn(0);
1177d763cef2SBarry Smith }
1178d763cef2SBarry Smith 
11794a2ae208SSatish Balay #undef __FUNCT__
11804a2ae208SSatish Balay #define __FUNCT__ "TSGetProblemType"
1181d763cef2SBarry Smith /*@C
1182d763cef2SBarry Smith    TSGetProblemType - Returns the problem type of a TS (timestepper) context.
1183d763cef2SBarry Smith 
1184d763cef2SBarry Smith    Not Collective
1185d763cef2SBarry Smith 
1186d763cef2SBarry Smith    Input Parameter:
1187d763cef2SBarry Smith .  ts   - The TS context obtained from TSCreate()
1188d763cef2SBarry Smith 
1189d763cef2SBarry Smith    Output Parameter:
1190d763cef2SBarry Smith .  type - The problem type, TS_LINEAR or TS_NONLINEAR
1191d763cef2SBarry Smith 
1192d763cef2SBarry Smith    Level: intermediate
1193d763cef2SBarry Smith 
1194d763cef2SBarry Smith    Contributed by: Matthew Knepley
1195d763cef2SBarry Smith 
1196d763cef2SBarry Smith .keywords: ts, get, type
1197d763cef2SBarry Smith 
1198d763cef2SBarry Smith @*/
1199d763cef2SBarry Smith int TSGetProblemType(TS ts,TSProblemType *type)
1200d763cef2SBarry Smith {
1201d763cef2SBarry Smith   PetscFunctionBegin;
1202d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1203d763cef2SBarry Smith   *type = ts->problem_type;
1204d763cef2SBarry Smith   PetscFunctionReturn(0);
1205d763cef2SBarry Smith }
1206d763cef2SBarry Smith 
12074a2ae208SSatish Balay #undef __FUNCT__
12084a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1209d763cef2SBarry Smith /*@C
1210d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1211d763cef2SBarry Smith    TS options in the database.
1212d763cef2SBarry Smith 
1213d763cef2SBarry Smith    Collective on TS
1214d763cef2SBarry Smith 
1215d763cef2SBarry Smith    Input Parameter:
1216d763cef2SBarry Smith +  ts     - The TS context
1217d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1218d763cef2SBarry Smith 
1219d763cef2SBarry Smith    Notes:
1220d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1221d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1222d763cef2SBarry Smith    hyphen.
1223d763cef2SBarry Smith 
1224d763cef2SBarry Smith    Contributed by: Matthew Knepley
1225d763cef2SBarry Smith 
1226d763cef2SBarry Smith    Level: advanced
1227d763cef2SBarry Smith 
1228d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1229d763cef2SBarry Smith 
1230d763cef2SBarry Smith .seealso: TSSetFromOptions()
1231d763cef2SBarry Smith 
1232d763cef2SBarry Smith @*/
1233d763cef2SBarry Smith int TSSetOptionsPrefix(TS ts,char *prefix)
1234d763cef2SBarry Smith {
1235d763cef2SBarry Smith   int ierr;
1236d763cef2SBarry Smith 
1237d763cef2SBarry Smith   PetscFunctionBegin;
1238d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1239d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1240d763cef2SBarry Smith   switch(ts->problem_type) {
1241d763cef2SBarry Smith     case TS_NONLINEAR:
1242d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1243d763cef2SBarry Smith       break;
1244d763cef2SBarry Smith     case TS_LINEAR:
1245d763cef2SBarry Smith       ierr = SLESSetOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr);
1246d763cef2SBarry Smith       break;
1247d763cef2SBarry Smith   }
1248d763cef2SBarry Smith   PetscFunctionReturn(0);
1249d763cef2SBarry Smith }
1250d763cef2SBarry Smith 
1251d763cef2SBarry Smith 
12524a2ae208SSatish Balay #undef __FUNCT__
12534a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1254d763cef2SBarry Smith /*@C
1255d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1256d763cef2SBarry Smith    TS options in the database.
1257d763cef2SBarry Smith 
1258d763cef2SBarry Smith    Collective on TS
1259d763cef2SBarry Smith 
1260d763cef2SBarry Smith    Input Parameter:
1261d763cef2SBarry Smith +  ts     - The TS context
1262d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1263d763cef2SBarry Smith 
1264d763cef2SBarry Smith    Notes:
1265d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1266d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1267d763cef2SBarry Smith    hyphen.
1268d763cef2SBarry Smith 
1269d763cef2SBarry Smith    Contributed by: Matthew Knepley
1270d763cef2SBarry Smith 
1271d763cef2SBarry Smith    Level: advanced
1272d763cef2SBarry Smith 
1273d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1274d763cef2SBarry Smith 
1275d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1276d763cef2SBarry Smith 
1277d763cef2SBarry Smith @*/
1278d763cef2SBarry Smith int TSAppendOptionsPrefix(TS ts,char *prefix)
1279d763cef2SBarry Smith {
1280d763cef2SBarry Smith   int ierr;
1281d763cef2SBarry Smith 
1282d763cef2SBarry Smith   PetscFunctionBegin;
1283d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1284d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1285d763cef2SBarry Smith   switch(ts->problem_type) {
1286d763cef2SBarry Smith     case TS_NONLINEAR:
1287d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1288d763cef2SBarry Smith       break;
1289d763cef2SBarry Smith     case TS_LINEAR:
1290d763cef2SBarry Smith       ierr = SLESAppendOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr);
1291d763cef2SBarry Smith       break;
1292d763cef2SBarry Smith   }
1293d763cef2SBarry Smith   PetscFunctionReturn(0);
1294d763cef2SBarry Smith }
1295d763cef2SBarry Smith 
12964a2ae208SSatish Balay #undef __FUNCT__
12974a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1298d763cef2SBarry Smith /*@C
1299d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1300d763cef2SBarry Smith    TS options in the database.
1301d763cef2SBarry Smith 
1302d763cef2SBarry Smith    Not Collective
1303d763cef2SBarry Smith 
1304d763cef2SBarry Smith    Input Parameter:
1305d763cef2SBarry Smith .  ts - The TS context
1306d763cef2SBarry Smith 
1307d763cef2SBarry Smith    Output Parameter:
1308d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1309d763cef2SBarry Smith 
1310d763cef2SBarry Smith    Contributed by: Matthew Knepley
1311d763cef2SBarry Smith 
1312d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1313d763cef2SBarry Smith    sufficient length to hold the prefix.
1314d763cef2SBarry Smith 
1315d763cef2SBarry Smith    Level: intermediate
1316d763cef2SBarry Smith 
1317d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1318d763cef2SBarry Smith 
1319d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1320d763cef2SBarry Smith @*/
1321d763cef2SBarry Smith int TSGetOptionsPrefix(TS ts,char **prefix)
1322d763cef2SBarry Smith {
1323d763cef2SBarry Smith   int ierr;
1324d763cef2SBarry Smith 
1325d763cef2SBarry Smith   PetscFunctionBegin;
1326d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1327d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1328d763cef2SBarry Smith   PetscFunctionReturn(0);
1329d763cef2SBarry Smith }
1330d763cef2SBarry Smith 
13314a2ae208SSatish Balay #undef __FUNCT__
13324a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1333d763cef2SBarry Smith /*@C
1334d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1335d763cef2SBarry Smith 
1336d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1337d763cef2SBarry Smith 
1338d763cef2SBarry Smith    Input Parameter:
1339d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1340d763cef2SBarry Smith 
1341d763cef2SBarry Smith    Output Parameters:
1342d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1343d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1344d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1345d763cef2SBarry Smith 
1346d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1347d763cef2SBarry Smith 
1348d763cef2SBarry Smith    Contributed by: Matthew Knepley
1349d763cef2SBarry Smith 
1350d763cef2SBarry Smith    Level: intermediate
1351d763cef2SBarry Smith 
1352d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1353d763cef2SBarry Smith 
1354d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1355d763cef2SBarry Smith 
1356d763cef2SBarry Smith @*/
1357d763cef2SBarry Smith int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1358d763cef2SBarry Smith {
1359d763cef2SBarry Smith   PetscFunctionBegin;
1360d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1361d763cef2SBarry Smith   if (A)   *A = ts->A;
1362d763cef2SBarry Smith   if (M)   *M = ts->B;
1363d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1364d763cef2SBarry Smith   PetscFunctionReturn(0);
1365d763cef2SBarry Smith }
1366d763cef2SBarry Smith 
13674a2ae208SSatish Balay #undef __FUNCT__
13684a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1369d763cef2SBarry Smith /*@C
1370d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1371d763cef2SBarry Smith 
1372d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1373d763cef2SBarry Smith 
1374d763cef2SBarry Smith    Input Parameter:
1375d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1376d763cef2SBarry Smith 
1377d763cef2SBarry Smith    Output Parameters:
1378d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1379d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1380d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1381d763cef2SBarry Smith 
1382d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1383d763cef2SBarry Smith 
1384d763cef2SBarry Smith    Contributed by: Matthew Knepley
1385d763cef2SBarry Smith 
1386d763cef2SBarry Smith    Level: intermediate
1387d763cef2SBarry Smith 
1388d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1389d763cef2SBarry Smith 
1390d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1391d763cef2SBarry Smith @*/
1392d763cef2SBarry Smith int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1393d763cef2SBarry Smith {
1394d763cef2SBarry Smith   int ierr;
1395d763cef2SBarry Smith 
1396d763cef2SBarry Smith   PetscFunctionBegin;
1397d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1398d763cef2SBarry Smith   PetscFunctionReturn(0);
1399d763cef2SBarry Smith }
1400d763cef2SBarry Smith 
1401d763cef2SBarry Smith /*MC
1402f1af5d2fSBarry Smith    TSRegisterDynamic - Adds a method to the timestepping solver package.
1403d763cef2SBarry Smith 
1404d763cef2SBarry Smith    Synopsis:
14053a7fca6bSBarry Smith    int TSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(TS))
1406d763cef2SBarry Smith 
1407d763cef2SBarry Smith    Not collective
1408d763cef2SBarry Smith 
1409d763cef2SBarry Smith    Input Parameters:
1410d763cef2SBarry Smith +  name_solver - name of a new user-defined solver
1411d763cef2SBarry Smith .  path - path (either absolute or relative) the library containing this solver
1412d763cef2SBarry Smith .  name_create - name of routine to create method context
1413d763cef2SBarry Smith -  routine_create - routine to create method context
1414d763cef2SBarry Smith 
1415d763cef2SBarry Smith    Notes:
1416f1af5d2fSBarry Smith    TSRegisterDynamic() may be called multiple times to add several user-defined solvers.
1417d763cef2SBarry Smith 
1418d763cef2SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1419d763cef2SBarry Smith    is ignored.
1420d763cef2SBarry Smith 
1421d763cef2SBarry Smith    Sample usage:
1422d763cef2SBarry Smith .vb
1423f1af5d2fSBarry Smith    TSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
1424d763cef2SBarry Smith               "MySolverCreate",MySolverCreate);
1425d763cef2SBarry Smith .ve
1426d763cef2SBarry Smith 
1427d763cef2SBarry Smith    Then, your solver can be chosen with the procedural interface via
1428d763cef2SBarry Smith $     TSSetType(ts,"my_solver")
1429d763cef2SBarry Smith    or at runtime via the option
1430d763cef2SBarry Smith $     -ts_type my_solver
1431d763cef2SBarry Smith 
1432d763cef2SBarry Smith    Level: advanced
1433d763cef2SBarry Smith 
14342aeb12f1SSatish Balay    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
14355ba88a07SLois Curfman McInnes    and others of the form ${any_environmental_variable} occuring in pathname will be
14365ba88a07SLois Curfman McInnes    replaced with appropriate values.
1437d763cef2SBarry Smith 
1438d763cef2SBarry Smith .keywords: TS, register
1439d763cef2SBarry Smith 
1440d763cef2SBarry Smith .seealso: TSRegisterAll(), TSRegisterDestroy()
1441d763cef2SBarry Smith M*/
1442d763cef2SBarry Smith 
14434a2ae208SSatish Balay #undef __FUNCT__
14444a2ae208SSatish Balay #define __FUNCT__ "TSRegister"
1445f1af5d2fSBarry Smith int TSRegister(char *sname,char *path,char *name,int (*function)(TS))
1446d763cef2SBarry Smith {
1447d763cef2SBarry Smith   char fullname[256];
1448549d3d68SSatish Balay   int  ierr;
1449d763cef2SBarry Smith 
1450d763cef2SBarry Smith   PetscFunctionBegin;
1451b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
14521713a123SBarry Smith   ierr = PetscFListAdd(&TSList,sname,fullname,(void (*)())function);CHKERRQ(ierr);
1453d763cef2SBarry Smith   PetscFunctionReturn(0);
1454d763cef2SBarry Smith }
14551713a123SBarry Smith 
14561713a123SBarry Smith #undef __FUNCT__
14571713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
14581713a123SBarry Smith /*@C
14591713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
14601713a123SBarry Smith    VecView() for the solution at each timestep
14611713a123SBarry Smith 
14621713a123SBarry Smith    Collective on TS
14631713a123SBarry Smith 
14641713a123SBarry Smith    Input Parameters:
14651713a123SBarry Smith +  ts - the TS context
14661713a123SBarry Smith .  step - current time-step
1467142b95e3SSatish Balay .  ptime - current time
14681713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
14691713a123SBarry Smith 
14701713a123SBarry Smith    Level: intermediate
14711713a123SBarry Smith 
14721713a123SBarry Smith .keywords: TS,  vector, monitor, view
14731713a123SBarry Smith 
14741713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
14751713a123SBarry Smith @*/
1476142b95e3SSatish Balay int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
14771713a123SBarry Smith {
14781713a123SBarry Smith   int         ierr;
14791713a123SBarry Smith   PetscViewer viewer = (PetscViewer) dummy;
14801713a123SBarry Smith 
14811713a123SBarry Smith   PetscFunctionBegin;
14821713a123SBarry Smith   if (!viewer) {
14831713a123SBarry Smith     MPI_Comm comm;
14841713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
14851713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
14861713a123SBarry Smith   }
14871713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
14881713a123SBarry Smith   PetscFunctionReturn(0);
14891713a123SBarry Smith }
14901713a123SBarry Smith 
14911713a123SBarry Smith 
14921713a123SBarry Smith 
1493