xref: /petsc/src/ts/interface/ts.c (revision d01c8aa3f8445be1931cee0defba57270449ff5a)
1 #define PETSCTS_DLL
2 
3 #include "include/private/tsimpl.h"        /*I "petscts.h"  I*/
4 
5 /* Logging support */
6 PetscCookie PETSCTS_DLLEXPORT TS_COOKIE = 0;
7 PetscEvent  TS_Step = 0, TS_PseudoComputeTimeStep = 0, TS_FunctionEval = 0, TS_JacobianEval = 0;
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "TSSetTypeFromOptions"
11 /*
12   TSSetTypeFromOptions - Sets the type of ts from user options.
13 
14   Collective on TS
15 
16   Input Parameter:
17 . ts - The ts
18 
19   Level: intermediate
20 
21 .keywords: TS, set, options, database, type
22 .seealso: TSSetFromOptions(), TSSetType()
23 */
24 static PetscErrorCode TSSetTypeFromOptions(TS ts)
25 {
26   PetscTruth     opt;
27   const char     *defaultType;
28   char           typeName[256];
29   PetscErrorCode ierr;
30 
31   PetscFunctionBegin;
32   if (ts->type_name) {
33     defaultType = ts->type_name;
34   } else {
35     defaultType = TS_EULER;
36   }
37 
38   if (!TSRegisterAllCalled) {
39     ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);
40   }
41   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
42   if (opt) {
43     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
44   } else {
45     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
46   }
47   PetscFunctionReturn(0);
48 }
49 
50 #undef __FUNCT__
51 #define __FUNCT__ "TSSetFromOptions"
52 /*@
53    TSSetFromOptions - Sets various TS parameters from user options.
54 
55    Collective on TS
56 
57    Input Parameter:
58 .  ts - the TS context obtained from TSCreate()
59 
60    Options Database Keys:
61 +  -ts_type <type> - TS_EULER, TS_BEULER, TS_SUNDIALS, TS_PSEUDO, TS_CRANK_NICHOLSON
62 .  -ts_max_steps maxsteps - maximum number of time-steps to take
63 .  -ts_max_time time - maximum time to compute to
64 .  -ts_dt dt - initial time step
65 .  -ts_monitor - print information at each timestep
66 -  -ts_monitor_draw - plot information at each timestep
67 
68    Level: beginner
69 
70 .keywords: TS, timestep, set, options, database
71 
72 .seealso: TSGetType
73 @*/
74 PetscErrorCode PETSCTS_DLLEXPORT TSSetFromOptions(TS ts)
75 {
76   PetscReal               dt;
77   PetscTruth              opt,flg;
78   PetscErrorCode          ierr;
79   PetscViewerASCIIMonitor monviewer;
80   char                    monfilename[PETSC_MAX_PATH_LEN];
81 
82   PetscFunctionBegin;
83   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
84   ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");CHKERRQ(ierr);
85 
86     /* Handle generic TS options */
87     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
88     ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
89     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
90     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
91     if (opt) {
92       ts->initial_time_step = ts->time_step = dt;
93     }
94 
95     /* Monitor options */
96     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
97     if (flg) {
98       ierr = PetscViewerASCIIMonitorCreate(ts->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
99       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
100     }
101     ierr = PetscOptionsName("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",&opt);CHKERRQ(ierr);
102     if (opt) {
103       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
104     }
105     ierr = PetscOptionsName("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",&opt);CHKERRQ(ierr);
106     if (opt) {
107       ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
108     }
109 
110     /* Handle TS type options */
111     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
112 
113     /* Handle specific TS options */
114     if (ts->ops->setfromoptions) {
115       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
116     }
117   ierr = PetscOptionsEnd();CHKERRQ(ierr);
118 
119   /* Handle subobject options */
120   switch(ts->problem_type) {
121     /* Should check for implicit/explicit */
122   case TS_LINEAR:
123     if (ts->ksp) {
124       ierr = KSPSetOperators(ts->ksp,ts->Arhs,ts->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
125       ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr);
126     }
127     break;
128   case TS_NONLINEAR:
129     if (ts->snes) {
130       /* this is a bit of a hack, but it gets the matrix information into SNES earlier
131          so that SNES and KSP have more information to pick reasonable defaults
132          before they allow users to set options */
133       ierr = SNESSetJacobian(ts->snes,ts->Arhs,ts->B,0,ts);CHKERRQ(ierr);
134       ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
135     }
136     break;
137   default:
138     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type);
139   }
140 
141   PetscFunctionReturn(0);
142 }
143 
144 #undef  __FUNCT__
145 #define __FUNCT__ "TSViewFromOptions"
146 /*@
147   TSViewFromOptions - This function visualizes the ts based upon user options.
148 
149   Collective on TS
150 
151   Input Parameter:
152 . ts - The ts
153 
154   Level: intermediate
155 
156 .keywords: TS, view, options, database
157 .seealso: TSSetFromOptions(), TSView()
158 @*/
159 PetscErrorCode PETSCTS_DLLEXPORT TSViewFromOptions(TS ts,const char title[])
160 {
161   PetscViewer    viewer;
162   PetscDraw      draw;
163   PetscTruth     opt;
164   char           fileName[PETSC_MAX_PATH_LEN];
165   PetscErrorCode ierr;
166 
167   PetscFunctionBegin;
168   ierr = PetscOptionsGetString(ts->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr);
169   if (opt && !PetscPreLoadingOn) {
170     ierr = PetscViewerASCIIOpen(ts->comm,fileName,&viewer);CHKERRQ(ierr);
171     ierr = TSView(ts, viewer);CHKERRQ(ierr);
172     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
173   }
174   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);CHKERRQ(ierr);
175   if (opt) {
176     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
177     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
178     if (title) {
179       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
180     } else {
181       ierr = PetscObjectName((PetscObject) ts);CHKERRQ(ierr);
182       ierr = PetscDrawSetTitle(draw, ts->name);CHKERRQ(ierr);
183     }
184     ierr = TSView(ts, viewer);CHKERRQ(ierr);
185     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
186     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
187     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
188   }
189   PetscFunctionReturn(0);
190 }
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "TSComputeRHSJacobian"
194 /*@
195    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
196       set with TSSetRHSJacobian().
197 
198    Collective on TS and Vec
199 
200    Input Parameters:
201 +  ts - the SNES context
202 .  t - current timestep
203 -  x - input vector
204 
205    Output Parameters:
206 +  A - Jacobian matrix
207 .  B - optional preconditioning matrix
208 -  flag - flag indicating matrix structure
209 
210    Notes:
211    Most users should not need to explicitly call this routine, as it
212    is used internally within the nonlinear solvers.
213 
214    See KSPSetOperators() for important information about setting the
215    flag parameter.
216 
217    TSComputeJacobian() is valid only for TS_NONLINEAR
218 
219    Level: developer
220 
221 .keywords: SNES, compute, Jacobian, matrix
222 
223 .seealso:  TSSetRHSJacobian(), KSPSetOperators()
224 @*/
225 PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
226 {
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
231   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
232   PetscCheckSameComm(ts,1,X,3);
233   if (ts->problem_type != TS_NONLINEAR) {
234     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
235   }
236   if (ts->ops->rhsjacobian) {
237     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
238     *flg = DIFFERENT_NONZERO_PATTERN;
239     PetscStackPush("TS user Jacobian function");
240     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
241     PetscStackPop;
242     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
243     /* make sure user returned a correct Jacobian and preconditioner */
244     PetscValidHeaderSpecific(*A,MAT_COOKIE,4);
245     PetscValidHeaderSpecific(*B,MAT_COOKIE,5);
246   } else {
247     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249     if (*A != *B) {
250       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252     }
253   }
254   PetscFunctionReturn(0);
255 }
256 
257 #undef __FUNCT__
258 #define __FUNCT__ "TSComputeRHSFunction"
259 /*
260    TSComputeRHSFunction - Evaluates the right-hand-side function.
261 
262    Note: If the user did not provide a function but merely a matrix,
263    this routine applies the matrix.
264 */
265 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
266 {
267   PetscErrorCode ierr;
268 
269   PetscFunctionBegin;
270   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
271   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
272   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
273 
274   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
275   if (ts->ops->rhsfunction) {
276     PetscStackPush("TS user right-hand-side function");
277     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
278     PetscStackPop;
279   } else {
280     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
281       MatStructure flg;
282       PetscStackPush("TS user right-hand-side matrix function");
283       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
284       PetscStackPop;
285     }
286     ierr = MatMult(ts->Arhs,x,y);CHKERRQ(ierr);
287   }
288 
289   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
290 
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "TSSetRHSFunction"
296 /*@C
297     TSSetRHSFunction - Sets the routine for evaluating the function,
298     F(t,u), where U_t = F(t,u).
299 
300     Collective on TS
301 
302     Input Parameters:
303 +   ts - the TS context obtained from TSCreate()
304 .   f - routine for evaluating the right-hand-side function
305 -   ctx - [optional] user-defined context for private data for the
306           function evaluation routine (may be PETSC_NULL)
307 
308     Calling sequence of func:
309 $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
310 
311 +   t - current timestep
312 .   u - input vector
313 .   F - function vector
314 -   ctx - [optional] user-defined function context
315 
316     Important:
317     The user MUST call either this routine or TSSetMatrices().
318 
319     Level: beginner
320 
321 .keywords: TS, timestep, set, right-hand-side, function
322 
323 .seealso: TSSetMatrices()
324 @*/
325 PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
326 {
327   PetscFunctionBegin;
328 
329   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
330   if (ts->problem_type == TS_LINEAR) {
331     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
332   }
333   ts->ops->rhsfunction = f;
334   ts->funP             = ctx;
335   PetscFunctionReturn(0);
336 }
337 
338 #undef __FUNCT__
339 #define __FUNCT__ "TSSetMatrices"
340 /*@C
341    TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs,
342    where Alhs(t) U_t = Arhs(t) U.
343 
344    Collective on TS
345 
346    Input Parameters:
347 +  ts   - the TS context obtained from TSCreate()
348 .  Arhs - matrix
349 .  frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
350           if Arhs is not a function of t.
351 .  Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix.
352 .  flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
353           if Alhs is not a function of t.
354 .  flag - flag indicating information about the matrix structure of Arhs and Alhs.
355           The available options are
356             SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs
357             DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs
358 -  ctx  - [optional] user-defined context for private data for the
359           matrix evaluation routine (may be PETSC_NULL)
360 
361    Calling sequence of func:
362 $     func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
363 
364 +  t - current timestep
365 .  A - matrix A, where U_t = A(t) U
366 .  B - preconditioner matrix, usually the same as A
367 .  flag - flag indicating information about the preconditioner matrix
368           structure (same as flag in KSPSetOperators())
369 -  ctx - [optional] user-defined context for matrix evaluation routine
370 
371    Notes:
372    The routine func() takes Mat* as the matrix arguments rather than Mat.
373    This allows the matrix evaluation routine to replace Arhs or Alhs with a
374    completely new new matrix structure (not just different matrix elements)
375    when appropriate, for instance, if the nonzero structure is changing
376    throughout the global iterations.
377 
378    Important:
379    The user MUST call either this routine or TSSetRHSFunction().
380 
381    Level: beginner
382 
383 .keywords: TS, timestep, set, matrix
384 
385 .seealso: TSSetRHSFunction()
386 @*/
387 PetscErrorCode PETSCTS_DLLEXPORT TSSetMatrices(TS ts,Mat Arhs,PetscErrorCode (*frhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat Alhs,PetscErrorCode (*flhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure flag,void *ctx)
388 {
389   PetscFunctionBegin;
390   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
391   PetscValidHeaderSpecific(Arhs,MAT_COOKIE,2);
392   PetscCheckSameComm(ts,1,Arhs,2);
393   if (Alhs) PetscCheckSameComm(ts,1,Alhs,4);
394   if (ts->problem_type == TS_NONLINEAR)
395     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
396 
397   ts->ops->rhsmatrix = frhs;
398   ts->ops->lhsmatrix = flhs;
399   ts->jacP           = ctx;
400   ts->Arhs           = Arhs;
401   ts->Alhs           = Alhs;
402   ts->matflg         = flag;
403   PetscFunctionReturn(0);
404 }
405 
406 #undef __FUNCT__
407 #define __FUNCT__ "TSGetMatrices"
408 /*@C
409    TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep,
410    where Alhs(t) U_t = Arhs(t) U.
411 
412    Not Collective, but parallel objects are returned if TS is parallel
413 
414    Input Parameter:
415 .  ts  - The TS context obtained from TSCreate()
416 
417    Output Parameters:
418 +  Arhs - The right-hand side matrix
419 .  Alhs - The left-hand side matrix
420 -  ctx - User-defined context for matrix evaluation routine
421 
422    Notes: You can pass in PETSC_NULL for any return argument you do not need.
423 
424    Level: intermediate
425 
426 .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
427 
428 .keywords: TS, timestep, get, matrix
429 
430 @*/
431 PetscErrorCode PETSCTS_DLLEXPORT TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx)
432 {
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
435   if (Arhs) *Arhs = ts->Arhs;
436   if (Alhs) *Alhs = ts->Alhs;
437   if (ctx)  *ctx = ts->jacP;
438   PetscFunctionReturn(0);
439 }
440 
441 #undef __FUNCT__
442 #define __FUNCT__ "TSSetRHSJacobian"
443 /*@C
444    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
445    where U_t = F(U,t), as well as the location to store the matrix.
446    Use TSSetMatrices() for linear problems.
447 
448    Collective on TS
449 
450    Input Parameters:
451 +  ts  - the TS context obtained from TSCreate()
452 .  A   - Jacobian matrix
453 .  B   - preconditioner matrix (usually same as A)
454 .  f   - the Jacobian evaluation routine
455 -  ctx - [optional] user-defined context for private data for the
456          Jacobian evaluation routine (may be PETSC_NULL)
457 
458    Calling sequence of func:
459 $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
460 
461 +  t - current timestep
462 .  u - input vector
463 .  A - matrix A, where U_t = A(t)u
464 .  B - preconditioner matrix, usually the same as A
465 .  flag - flag indicating information about the preconditioner matrix
466           structure (same as flag in KSPSetOperators())
467 -  ctx - [optional] user-defined context for matrix evaluation routine
468 
469    Notes:
470    See KSPSetOperators() for important information about setting the flag
471    output parameter in the routine func().  Be sure to read this information!
472 
473    The routine func() takes Mat * as the matrix arguments rather than Mat.
474    This allows the matrix evaluation routine to replace A and/or B with a
475    completely new new matrix structure (not just different matrix elements)
476    when appropriate, for instance, if the nonzero structure is changing
477    throughout the global iterations.
478 
479    Level: beginner
480 
481 .keywords: TS, timestep, set, right-hand-side, Jacobian
482 
483 .seealso: TSDefaultComputeJacobianColor(),
484           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices()
485 
486 @*/
487 PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
488 {
489   PetscFunctionBegin;
490   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
491   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
492   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
493   PetscCheckSameComm(ts,1,A,2);
494   PetscCheckSameComm(ts,1,B,3);
495   if (ts->problem_type != TS_NONLINEAR) {
496     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()");
497   }
498 
499   ts->ops->rhsjacobian = f;
500   ts->jacP             = ctx;
501   ts->Arhs             = A;
502   ts->B                = B;
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "TSView"
508 /*@C
509     TSView - Prints the TS data structure.
510 
511     Collective on TS
512 
513     Input Parameters:
514 +   ts - the TS context obtained from TSCreate()
515 -   viewer - visualization context
516 
517     Options Database Key:
518 .   -ts_view - calls TSView() at end of TSStep()
519 
520     Notes:
521     The available visualization contexts include
522 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
523 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
524          output where only the first processor opens
525          the file.  All other processors send their
526          data to the first processor to print.
527 
528     The user can open an alternative visualization context with
529     PetscViewerASCIIOpen() - output to a specified file.
530 
531     Level: beginner
532 
533 .keywords: TS, timestep, view
534 
535 .seealso: PetscViewerASCIIOpen()
536 @*/
537 PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer)
538 {
539   PetscErrorCode ierr;
540   char           *type;
541   PetscTruth     iascii,isstring;
542 
543   PetscFunctionBegin;
544   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
545   if (!viewer) {
546     ierr = PetscViewerASCIIGetStdout(ts->comm,&viewer);CHKERRQ(ierr);
547   }
548   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
549   PetscCheckSameComm(ts,1,viewer,2);
550 
551   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
552   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
553   if (iascii) {
554     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
555     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
556     if (type) {
557       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
558     } else {
559       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
560     }
561     if (ts->ops->view) {
562       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
563       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
564       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
565     }
566     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
567     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
568     if (ts->problem_type == TS_NONLINEAR) {
569       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
570     }
571     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
572   } else if (isstring) {
573     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
574     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
575   }
576   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
577   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
578   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
579   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
580   PetscFunctionReturn(0);
581 }
582 
583 
584 #undef __FUNCT__
585 #define __FUNCT__ "TSSetApplicationContext"
586 /*@C
587    TSSetApplicationContext - Sets an optional user-defined context for
588    the timesteppers.
589 
590    Collective on TS
591 
592    Input Parameters:
593 +  ts - the TS context obtained from TSCreate()
594 -  usrP - optional user context
595 
596    Level: intermediate
597 
598 .keywords: TS, timestep, set, application, context
599 
600 .seealso: TSGetApplicationContext()
601 @*/
602 PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP)
603 {
604   PetscFunctionBegin;
605   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
606   ts->user = usrP;
607   PetscFunctionReturn(0);
608 }
609 
610 #undef __FUNCT__
611 #define __FUNCT__ "TSGetApplicationContext"
612 /*@C
613     TSGetApplicationContext - Gets the user-defined context for the
614     timestepper.
615 
616     Not Collective
617 
618     Input Parameter:
619 .   ts - the TS context obtained from TSCreate()
620 
621     Output Parameter:
622 .   usrP - user context
623 
624     Level: intermediate
625 
626 .keywords: TS, timestep, get, application, context
627 
628 .seealso: TSSetApplicationContext()
629 @*/
630 PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP)
631 {
632   PetscFunctionBegin;
633   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
634   *usrP = ts->user;
635   PetscFunctionReturn(0);
636 }
637 
638 #undef __FUNCT__
639 #define __FUNCT__ "TSGetTimeStepNumber"
640 /*@
641    TSGetTimeStepNumber - Gets the current number of timesteps.
642 
643    Not Collective
644 
645    Input Parameter:
646 .  ts - the TS context obtained from TSCreate()
647 
648    Output Parameter:
649 .  iter - number steps so far
650 
651    Level: intermediate
652 
653 .keywords: TS, timestep, get, iteration, number
654 @*/
655 PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter)
656 {
657   PetscFunctionBegin;
658   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
659   PetscValidIntPointer(iter,2);
660   *iter = ts->steps;
661   PetscFunctionReturn(0);
662 }
663 
664 #undef __FUNCT__
665 #define __FUNCT__ "TSSetInitialTimeStep"
666 /*@
667    TSSetInitialTimeStep - Sets the initial timestep to be used,
668    as well as the initial time.
669 
670    Collective on TS
671 
672    Input Parameters:
673 +  ts - the TS context obtained from TSCreate()
674 .  initial_time - the initial time
675 -  time_step - the size of the timestep
676 
677    Level: intermediate
678 
679 .seealso: TSSetTimeStep(), TSGetTimeStep()
680 
681 .keywords: TS, set, initial, timestep
682 @*/
683 PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
684 {
685   PetscFunctionBegin;
686   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
687   ts->time_step         = time_step;
688   ts->initial_time_step = time_step;
689   ts->ptime             = initial_time;
690   PetscFunctionReturn(0);
691 }
692 
693 #undef __FUNCT__
694 #define __FUNCT__ "TSSetTimeStep"
695 /*@
696    TSSetTimeStep - Allows one to reset the timestep at any time,
697    useful for simple pseudo-timestepping codes.
698 
699    Collective on TS
700 
701    Input Parameters:
702 +  ts - the TS context obtained from TSCreate()
703 -  time_step - the size of the timestep
704 
705    Level: intermediate
706 
707 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
708 
709 .keywords: TS, set, timestep
710 @*/
711 PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step)
712 {
713   PetscFunctionBegin;
714   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
715   ts->time_step = time_step;
716   PetscFunctionReturn(0);
717 }
718 
719 #undef __FUNCT__
720 #define __FUNCT__ "TSGetTimeStep"
721 /*@
722    TSGetTimeStep - Gets the current timestep size.
723 
724    Not Collective
725 
726    Input Parameter:
727 .  ts - the TS context obtained from TSCreate()
728 
729    Output Parameter:
730 .  dt - the current timestep size
731 
732    Level: intermediate
733 
734 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
735 
736 .keywords: TS, get, timestep
737 @*/
738 PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt)
739 {
740   PetscFunctionBegin;
741   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
742   PetscValidDoublePointer(dt,2);
743   *dt = ts->time_step;
744   PetscFunctionReturn(0);
745 }
746 
747 #undef __FUNCT__
748 #define __FUNCT__ "TSGetSolution"
749 /*@
750    TSGetSolution - Returns the solution at the present timestep. It
751    is valid to call this routine inside the function that you are evaluating
752    in order to move to the new timestep. This vector not changed until
753    the solution at the next timestep has been calculated.
754 
755    Not Collective, but Vec returned is parallel if TS is parallel
756 
757    Input Parameter:
758 .  ts - the TS context obtained from TSCreate()
759 
760    Output Parameter:
761 .  v - the vector containing the solution
762 
763    Level: intermediate
764 
765 .seealso: TSGetTimeStep()
766 
767 .keywords: TS, timestep, get, solution
768 @*/
769 PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v)
770 {
771   PetscFunctionBegin;
772   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
773   PetscValidPointer(v,2);
774   *v = ts->vec_sol_always;
775   PetscFunctionReturn(0);
776 }
777 
778 /* ----- Routines to initialize and destroy a timestepper ---- */
779 #undef __FUNCT__
780 #define __FUNCT__ "TSSetProblemType"
781 /*@
782   TSSetProblemType - Sets the type of problem to be solved.
783 
784   Not collective
785 
786   Input Parameters:
787 + ts   - The TS
788 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
789 .vb
790          U_t = A U
791          U_t = A(t) U
792          U_t = F(t,U)
793 .ve
794 
795    Level: beginner
796 
797 .keywords: TS, problem type
798 .seealso: TSSetUp(), TSProblemType, TS
799 @*/
800 PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type)
801 {
802   PetscFunctionBegin;
803   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
804   ts->problem_type = type;
805   PetscFunctionReturn(0);
806 }
807 
808 #undef __FUNCT__
809 #define __FUNCT__ "TSGetProblemType"
810 /*@C
811   TSGetProblemType - Gets the type of problem to be solved.
812 
813   Not collective
814 
815   Input Parameter:
816 . ts   - The TS
817 
818   Output Parameter:
819 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
820 .vb
821          U_t = A U
822          U_t = A(t) U
823          U_t = F(t,U)
824 .ve
825 
826    Level: beginner
827 
828 .keywords: TS, problem type
829 .seealso: TSSetUp(), TSProblemType, TS
830 @*/
831 PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type)
832 {
833   PetscFunctionBegin;
834   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
835   PetscValidIntPointer(type,2);
836   *type = ts->problem_type;
837   PetscFunctionReturn(0);
838 }
839 
840 #undef __FUNCT__
841 #define __FUNCT__ "TSSetUp"
842 /*@
843    TSSetUp - Sets up the internal data structures for the later use
844    of a timestepper.
845 
846    Collective on TS
847 
848    Input Parameter:
849 .  ts - the TS context obtained from TSCreate()
850 
851    Notes:
852    For basic use of the TS solvers the user need not explicitly call
853    TSSetUp(), since these actions will automatically occur during
854    the call to TSStep().  However, if one wishes to control this
855    phase separately, TSSetUp() should be called after TSCreate()
856    and optional routines of the form TSSetXXX(), but before TSStep().
857 
858    Level: advanced
859 
860 .keywords: TS, timestep, setup
861 
862 .seealso: TSCreate(), TSStep(), TSDestroy()
863 @*/
864 PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts)
865 {
866   PetscErrorCode ierr;
867 
868   PetscFunctionBegin;
869   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
870   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
871   if (!ts->type_name) {
872     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
873   }
874   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
875   ts->setupcalled = 1;
876   PetscFunctionReturn(0);
877 }
878 
879 #undef __FUNCT__
880 #define __FUNCT__ "TSDestroy"
881 /*@
882    TSDestroy - Destroys the timestepper context that was created
883    with TSCreate().
884 
885    Collective on TS
886 
887    Input Parameter:
888 .  ts - the TS context obtained from TSCreate()
889 
890    Level: beginner
891 
892 .keywords: TS, timestepper, destroy
893 
894 .seealso: TSCreate(), TSSetUp(), TSSolve()
895 @*/
896 PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts)
897 {
898   PetscErrorCode ierr;
899 
900   PetscFunctionBegin;
901   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
902   if (--ts->refct > 0) PetscFunctionReturn(0);
903 
904   /* if memory was published with AMS then destroy it */
905   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
906   if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr)}
907   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
908   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
909   if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);}
910   ierr = TSMonitorCancel(ts);CHKERRQ(ierr);
911   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
912   PetscFunctionReturn(0);
913 }
914 
915 #undef __FUNCT__
916 #define __FUNCT__ "TSGetSNES"
917 /*@
918    TSGetSNES - Returns the SNES (nonlinear solver) associated with
919    a TS (timestepper) context. Valid only for nonlinear problems.
920 
921    Not Collective, but SNES is parallel if TS is parallel
922 
923    Input Parameter:
924 .  ts - the TS context obtained from TSCreate()
925 
926    Output Parameter:
927 .  snes - the nonlinear solver context
928 
929    Notes:
930    The user can then directly manipulate the SNES context to set various
931    options, etc.  Likewise, the user can then extract and manipulate the
932    KSP, KSP, and PC contexts as well.
933 
934    TSGetSNES() does not work for integrators that do not use SNES; in
935    this case TSGetSNES() returns PETSC_NULL in snes.
936 
937    Level: beginner
938 
939 .keywords: timestep, get, SNES
940 @*/
941 PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes)
942 {
943   PetscFunctionBegin;
944   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
945   PetscValidPointer(snes,2);
946   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
947   *snes = ts->snes;
948   PetscFunctionReturn(0);
949 }
950 
951 #undef __FUNCT__
952 #define __FUNCT__ "TSGetKSP"
953 /*@
954    TSGetKSP - Returns the KSP (linear solver) associated with
955    a TS (timestepper) context.
956 
957    Not Collective, but KSP is parallel if TS is parallel
958 
959    Input Parameter:
960 .  ts - the TS context obtained from TSCreate()
961 
962    Output Parameter:
963 .  ksp - the nonlinear solver context
964 
965    Notes:
966    The user can then directly manipulate the KSP context to set various
967    options, etc.  Likewise, the user can then extract and manipulate the
968    KSP and PC contexts as well.
969 
970    TSGetKSP() does not work for integrators that do not use KSP;
971    in this case TSGetKSP() returns PETSC_NULL in ksp.
972 
973    Level: beginner
974 
975 .keywords: timestep, get, KSP
976 @*/
977 PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp)
978 {
979   PetscFunctionBegin;
980   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
981   PetscValidPointer(ksp,2);
982   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
983   *ksp = ts->ksp;
984   PetscFunctionReturn(0);
985 }
986 
987 /* ----------- Routines to set solver parameters ---------- */
988 
989 #undef __FUNCT__
990 #define __FUNCT__ "TSGetDuration"
991 /*@
992    TSGetDuration - Gets the maximum number of timesteps to use and
993    maximum time for iteration.
994 
995    Collective on TS
996 
997    Input Parameters:
998 +  ts       - the TS context obtained from TSCreate()
999 .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1000 -  maxtime  - final time to iterate to, or PETSC_NULL
1001 
1002    Level: intermediate
1003 
1004 .keywords: TS, timestep, get, maximum, iterations, time
1005 @*/
1006 PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1007 {
1008   PetscFunctionBegin;
1009   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1010   if (maxsteps) {
1011     PetscValidIntPointer(maxsteps,2);
1012     *maxsteps = ts->max_steps;
1013   }
1014   if (maxtime ) {
1015     PetscValidScalarPointer(maxtime,3);
1016     *maxtime  = ts->max_time;
1017   }
1018   PetscFunctionReturn(0);
1019 }
1020 
1021 #undef __FUNCT__
1022 #define __FUNCT__ "TSSetDuration"
1023 /*@
1024    TSSetDuration - Sets the maximum number of timesteps to use and
1025    maximum time for iteration.
1026 
1027    Collective on TS
1028 
1029    Input Parameters:
1030 +  ts - the TS context obtained from TSCreate()
1031 .  maxsteps - maximum number of iterations to use
1032 -  maxtime - final time to iterate to
1033 
1034    Options Database Keys:
1035 .  -ts_max_steps <maxsteps> - Sets maxsteps
1036 .  -ts_max_time <maxtime> - Sets maxtime
1037 
1038    Notes:
1039    The default maximum number of iterations is 5000. Default time is 5.0
1040 
1041    Level: intermediate
1042 
1043 .keywords: TS, timestep, set, maximum, iterations
1044 @*/
1045 PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1046 {
1047   PetscFunctionBegin;
1048   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1049   ts->max_steps = maxsteps;
1050   ts->max_time  = maxtime;
1051   PetscFunctionReturn(0);
1052 }
1053 
1054 #undef __FUNCT__
1055 #define __FUNCT__ "TSSetSolution"
1056 /*@
1057    TSSetSolution - Sets the initial solution vector
1058    for use by the TS routines.
1059 
1060    Collective on TS and Vec
1061 
1062    Input Parameters:
1063 +  ts - the TS context obtained from TSCreate()
1064 -  x - the solution vector
1065 
1066    Level: beginner
1067 
1068 .keywords: TS, timestep, set, solution, initial conditions
1069 @*/
1070 PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x)
1071 {
1072   PetscFunctionBegin;
1073   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1074   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1075   ts->vec_sol        = ts->vec_sol_always = x;
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 #undef __FUNCT__
1080 #define __FUNCT__ "TSSetPreStep"
1081 /*@C
1082   TSSetPreStep - Sets the general-purpose function
1083   called once at the beginning of time stepping.
1084 
1085   Collective on TS
1086 
1087   Input Parameters:
1088 + ts   - The TS context obtained from TSCreate()
1089 - func - The function
1090 
1091   Calling sequence of func:
1092 . func (TS ts);
1093 
1094   Level: intermediate
1095 
1096 .keywords: TS, timestep
1097 @*/
1098 PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1099 {
1100   PetscFunctionBegin;
1101   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1102   ts->ops->prestep = func;
1103   PetscFunctionReturn(0);
1104 }
1105 
1106 #undef __FUNCT__
1107 #define __FUNCT__ "TSDefaultPreStep"
1108 /*@
1109   TSDefaultPreStep - The default pre-stepping function which does nothing.
1110 
1111   Collective on TS
1112 
1113   Input Parameters:
1114 . ts  - The TS context obtained from TSCreate()
1115 
1116   Level: developer
1117 
1118 .keywords: TS, timestep
1119 @*/
1120 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts)
1121 {
1122   PetscFunctionBegin;
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 #undef __FUNCT__
1127 #define __FUNCT__ "TSSetUpdate"
1128 /*@C
1129   TSSetUpdate - Sets the general-purpose update function called
1130   at the beginning of every time step. This function can change
1131   the time step.
1132 
1133   Collective on TS
1134 
1135   Input Parameters:
1136 + ts   - The TS context obtained from TSCreate()
1137 - func - The function
1138 
1139   Calling sequence of func:
1140 . func (TS ts, double t, double *dt);
1141 
1142 + t   - The current time
1143 - dt  - The current time step
1144 
1145   Level: intermediate
1146 
1147 .keywords: TS, update, timestep
1148 @*/
1149 PetscErrorCode PETSCTS_DLLEXPORT TSSetUpdate(TS ts, PetscErrorCode (*func)(TS, PetscReal, PetscReal *))
1150 {
1151   PetscFunctionBegin;
1152   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1153   ts->ops->update = func;
1154   PetscFunctionReturn(0);
1155 }
1156 
1157 #undef __FUNCT__
1158 #define __FUNCT__ "TSDefaultUpdate"
1159 /*@
1160   TSDefaultUpdate - The default update function which does nothing.
1161 
1162   Collective on TS
1163 
1164   Input Parameters:
1165 + ts  - The TS context obtained from TSCreate()
1166 - t   - The current time
1167 
1168   Output Parameters:
1169 . dt  - The current time step
1170 
1171   Level: developer
1172 
1173 .keywords: TS, update, timestep
1174 @*/
1175 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1176 {
1177   PetscFunctionBegin;
1178   PetscFunctionReturn(0);
1179 }
1180 
1181 #undef __FUNCT__
1182 #define __FUNCT__ "TSSetPostStep"
1183 /*@C
1184   TSSetPostStep - Sets the general-purpose function
1185   called once at the end of time stepping.
1186 
1187   Collective on TS
1188 
1189   Input Parameters:
1190 + ts   - The TS context obtained from TSCreate()
1191 - func - The function
1192 
1193   Calling sequence of func:
1194 . func (TS ts);
1195 
1196   Level: intermediate
1197 
1198 .keywords: TS, timestep
1199 @*/
1200 PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1201 {
1202   PetscFunctionBegin;
1203   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1204   ts->ops->poststep = func;
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "TSDefaultPostStep"
1210 /*@
1211   TSDefaultPostStep - The default post-stepping function which does nothing.
1212 
1213   Collective on TS
1214 
1215   Input Parameters:
1216 . ts  - The TS context obtained from TSCreate()
1217 
1218   Level: developer
1219 
1220 .keywords: TS, timestep
1221 @*/
1222 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts)
1223 {
1224   PetscFunctionBegin;
1225   PetscFunctionReturn(0);
1226 }
1227 
1228 /* ------------ Routines to set performance monitoring options ----------- */
1229 
1230 #undef __FUNCT__
1231 #define __FUNCT__ "TSMonitorSet"
1232 /*@C
1233    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1234    timestep to display the iteration's  progress.
1235 
1236    Collective on TS
1237 
1238    Input Parameters:
1239 +  ts - the TS context obtained from TSCreate()
1240 .  func - monitoring routine
1241 .  mctx - [optional] user-defined context for private data for the
1242              monitor routine (use PETSC_NULL if no context is desired)
1243 -  monitordestroy - [optional] routine that frees monitor context
1244           (may be PETSC_NULL)
1245 
1246    Calling sequence of func:
1247 $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1248 
1249 +    ts - the TS context
1250 .    steps - iteration number
1251 .    time - current time
1252 .    x - current iterate
1253 -    mctx - [optional] monitoring context
1254 
1255    Notes:
1256    This routine adds an additional monitor to the list of monitors that
1257    already has been loaded.
1258 
1259    Level: intermediate
1260 
1261 .keywords: TS, timestep, set, monitor
1262 
1263 .seealso: TSMonitorDefault(), TSMonitorCancel()
1264 @*/
1265 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1266 {
1267   PetscFunctionBegin;
1268   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1269   if (ts->numbermonitors >= MAXTSMONITORS) {
1270     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1271   }
1272   ts->monitor[ts->numbermonitors]           = monitor;
1273   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1274   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 #undef __FUNCT__
1279 #define __FUNCT__ "TSMonitorCancel"
1280 /*@C
1281    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1282 
1283    Collective on TS
1284 
1285    Input Parameters:
1286 .  ts - the TS context obtained from TSCreate()
1287 
1288    Notes:
1289    There is no way to remove a single, specific monitor.
1290 
1291    Level: intermediate
1292 
1293 .keywords: TS, timestep, set, monitor
1294 
1295 .seealso: TSMonitorDefault(), TSMonitorSet()
1296 @*/
1297 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts)
1298 {
1299   PetscErrorCode ierr;
1300   PetscInt       i;
1301 
1302   PetscFunctionBegin;
1303   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1304   for (i=0; i<ts->numbermonitors; i++) {
1305     if (ts->mdestroy[i]) {
1306       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
1307     }
1308   }
1309   ts->numbermonitors = 0;
1310   PetscFunctionReturn(0);
1311 }
1312 
1313 #undef __FUNCT__
1314 #define __FUNCT__ "TSMonitorDefault"
1315 /*@
1316    TSMonitorDefault - Sets the Default monitor
1317 @*/
1318 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1319 {
1320   PetscErrorCode          ierr;
1321   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
1322 
1323   PetscFunctionBegin;
1324   if (!ctx) {
1325     ierr = PetscViewerASCIIMonitorCreate(ts->comm,"stdout",0,&viewer);CHKERRQ(ierr);
1326   }
1327   ierr = PetscViewerASCIIMonitorPrintf(viewer,"timestep %D dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1328   if (!ctx) {
1329     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
1330   }
1331   PetscFunctionReturn(0);
1332 }
1333 
1334 #undef __FUNCT__
1335 #define __FUNCT__ "TSStep"
1336 /*@
1337    TSStep - Steps the requested number of timesteps.
1338 
1339    Collective on TS
1340 
1341    Input Parameter:
1342 .  ts - the TS context obtained from TSCreate()
1343 
1344    Output Parameters:
1345 +  steps - number of iterations until termination
1346 -  ptime - time until termination
1347 
1348    Level: beginner
1349 
1350 .keywords: TS, timestep, solve
1351 
1352 .seealso: TSCreate(), TSSetUp(), TSDestroy()
1353 @*/
1354 PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1355 {
1356   PetscErrorCode ierr;
1357 
1358   PetscFunctionBegin;
1359   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1360   if (!ts->setupcalled) {
1361     ierr = TSSetUp(ts);CHKERRQ(ierr);
1362   }
1363 
1364   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1365   ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1366   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1367   ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1368   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1369 
1370   if (!PetscPreLoadingOn) {
1371     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1372   }
1373   PetscFunctionReturn(0);
1374 }
1375 
1376 #undef __FUNCT__
1377 #define __FUNCT__ "TSSolve"
1378 /*@
1379    TSSolve - Steps the requested number of timesteps.
1380 
1381    Collective on TS
1382 
1383    Input Parameter:
1384 +  ts - the TS context obtained from TSCreate()
1385 -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
1386 
1387    Level: beginner
1388 
1389 .keywords: TS, timestep, solve
1390 
1391 .seealso: TSCreate(), TSSetSolution(), TSStep()
1392 @*/
1393 PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x)
1394 {
1395   PetscInt       steps;
1396   PetscReal      ptime;
1397   PetscErrorCode ierr;
1398   PetscFunctionBegin;
1399   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1400   /* set solution vector if provided */
1401   if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); }
1402   /* reset time step and iteration counters */
1403   ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0;
1404   /* steps the requested number of timesteps. */
1405   ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr);
1406   PetscFunctionReturn(0);
1407 }
1408 
1409 #undef __FUNCT__
1410 #define __FUNCT__ "TSMonitor"
1411 /*
1412      Runs the user provided monitor routines, if they exists.
1413 */
1414 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1415 {
1416   PetscErrorCode ierr;
1417   PetscInt i,n = ts->numbermonitors;
1418 
1419   PetscFunctionBegin;
1420   for (i=0; i<n; i++) {
1421     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1422   }
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 /* ------------------------------------------------------------------------*/
1427 
1428 #undef __FUNCT__
1429 #define __FUNCT__ "TSMonitorLGCreate"
1430 /*@C
1431    TSMonitorLGCreate - Creates a line graph context for use with
1432    TS to monitor convergence of preconditioned residual norms.
1433 
1434    Collective on TS
1435 
1436    Input Parameters:
1437 +  host - the X display to open, or null for the local machine
1438 .  label - the title to put in the title bar
1439 .  x, y - the screen coordinates of the upper left coordinate of the window
1440 -  m, n - the screen width and height in pixels
1441 
1442    Output Parameter:
1443 .  draw - the drawing context
1444 
1445    Options Database Key:
1446 .  -ts_monitor_draw - automatically sets line graph monitor
1447 
1448    Notes:
1449    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1450 
1451    Level: intermediate
1452 
1453 .keywords: TS, monitor, line graph, residual, seealso
1454 
1455 .seealso: TSMonitorLGDestroy(), TSMonitorSet()
1456 
1457 @*/
1458 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1459 {
1460   PetscDraw      win;
1461   PetscErrorCode ierr;
1462 
1463   PetscFunctionBegin;
1464   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1465   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1466   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1467   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1468 
1469   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1470   PetscFunctionReturn(0);
1471 }
1472 
1473 #undef __FUNCT__
1474 #define __FUNCT__ "TSMonitorLG"
1475 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1476 {
1477   PetscDrawLG    lg = (PetscDrawLG) monctx;
1478   PetscReal      x,y = ptime;
1479   PetscErrorCode ierr;
1480 
1481   PetscFunctionBegin;
1482   if (!monctx) {
1483     MPI_Comm    comm;
1484     PetscViewer viewer;
1485 
1486     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1487     viewer = PETSC_VIEWER_DRAW_(comm);
1488     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
1489   }
1490 
1491   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1492   x = (PetscReal)n;
1493   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1494   if (n < 20 || (n % 5)) {
1495     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1496   }
1497   PetscFunctionReturn(0);
1498 }
1499 
1500 #undef __FUNCT__
1501 #define __FUNCT__ "TSMonitorLGDestroy"
1502 /*@C
1503    TSMonitorLGDestroy - Destroys a line graph context that was created
1504    with TSMonitorLGCreate().
1505 
1506    Collective on PetscDrawLG
1507 
1508    Input Parameter:
1509 .  draw - the drawing context
1510 
1511    Level: intermediate
1512 
1513 .keywords: TS, monitor, line graph, destroy
1514 
1515 .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1516 @*/
1517 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg)
1518 {
1519   PetscDraw      draw;
1520   PetscErrorCode ierr;
1521 
1522   PetscFunctionBegin;
1523   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1524   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1525   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1526   PetscFunctionReturn(0);
1527 }
1528 
1529 #undef __FUNCT__
1530 #define __FUNCT__ "TSGetTime"
1531 /*@
1532    TSGetTime - Gets the current time.
1533 
1534    Not Collective
1535 
1536    Input Parameter:
1537 .  ts - the TS context obtained from TSCreate()
1538 
1539    Output Parameter:
1540 .  t  - the current time
1541 
1542    Contributed by: Matthew Knepley
1543 
1544    Level: beginner
1545 
1546 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1547 
1548 .keywords: TS, get, time
1549 @*/
1550 PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t)
1551 {
1552   PetscFunctionBegin;
1553   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1554   PetscValidDoublePointer(t,2);
1555   *t = ts->ptime;
1556   PetscFunctionReturn(0);
1557 }
1558 
1559 #undef __FUNCT__
1560 #define __FUNCT__ "TSSetTime"
1561 /*@
1562    TSSetTime - Allows one to reset the time.
1563 
1564    Collective on TS
1565 
1566    Input Parameters:
1567 +  ts - the TS context obtained from TSCreate()
1568 -  time - the time
1569 
1570    Level: intermediate
1571 
1572 .seealso: TSGetTime(), TSSetDuration()
1573 
1574 .keywords: TS, set, time
1575 @*/
1576 PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t)
1577 {
1578   PetscFunctionBegin;
1579   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1580   ts->ptime = t;
1581   PetscFunctionReturn(0);
1582 }
1583 
1584 #undef __FUNCT__
1585 #define __FUNCT__ "TSSetOptionsPrefix"
1586 /*@C
1587    TSSetOptionsPrefix - Sets the prefix used for searching for all
1588    TS options in the database.
1589 
1590    Collective on TS
1591 
1592    Input Parameter:
1593 +  ts     - The TS context
1594 -  prefix - The prefix to prepend to all option names
1595 
1596    Notes:
1597    A hyphen (-) must NOT be given at the beginning of the prefix name.
1598    The first character of all runtime options is AUTOMATICALLY the
1599    hyphen.
1600 
1601    Contributed by: Matthew Knepley
1602 
1603    Level: advanced
1604 
1605 .keywords: TS, set, options, prefix, database
1606 
1607 .seealso: TSSetFromOptions()
1608 
1609 @*/
1610 PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[])
1611 {
1612   PetscErrorCode ierr;
1613 
1614   PetscFunctionBegin;
1615   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1616   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1617   switch(ts->problem_type) {
1618     case TS_NONLINEAR:
1619       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1620       break;
1621     case TS_LINEAR:
1622       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1623       break;
1624   }
1625   PetscFunctionReturn(0);
1626 }
1627 
1628 
1629 #undef __FUNCT__
1630 #define __FUNCT__ "TSAppendOptionsPrefix"
1631 /*@C
1632    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1633    TS options in the database.
1634 
1635    Collective on TS
1636 
1637    Input Parameter:
1638 +  ts     - The TS context
1639 -  prefix - The prefix to prepend to all option names
1640 
1641    Notes:
1642    A hyphen (-) must NOT be given at the beginning of the prefix name.
1643    The first character of all runtime options is AUTOMATICALLY the
1644    hyphen.
1645 
1646    Contributed by: Matthew Knepley
1647 
1648    Level: advanced
1649 
1650 .keywords: TS, append, options, prefix, database
1651 
1652 .seealso: TSGetOptionsPrefix()
1653 
1654 @*/
1655 PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[])
1656 {
1657   PetscErrorCode ierr;
1658 
1659   PetscFunctionBegin;
1660   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1661   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1662   switch(ts->problem_type) {
1663     case TS_NONLINEAR:
1664       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1665       break;
1666     case TS_LINEAR:
1667       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1668       break;
1669   }
1670   PetscFunctionReturn(0);
1671 }
1672 
1673 #undef __FUNCT__
1674 #define __FUNCT__ "TSGetOptionsPrefix"
1675 /*@C
1676    TSGetOptionsPrefix - Sets the prefix used for searching for all
1677    TS options in the database.
1678 
1679    Not Collective
1680 
1681    Input Parameter:
1682 .  ts - The TS context
1683 
1684    Output Parameter:
1685 .  prefix - A pointer to the prefix string used
1686 
1687    Contributed by: Matthew Knepley
1688 
1689    Notes: On the fortran side, the user should pass in a string 'prifix' of
1690    sufficient length to hold the prefix.
1691 
1692    Level: intermediate
1693 
1694 .keywords: TS, get, options, prefix, database
1695 
1696 .seealso: TSAppendOptionsPrefix()
1697 @*/
1698 PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[])
1699 {
1700   PetscErrorCode ierr;
1701 
1702   PetscFunctionBegin;
1703   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1704   PetscValidPointer(prefix,2);
1705   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1706   PetscFunctionReturn(0);
1707 }
1708 
1709 #undef __FUNCT__
1710 #define __FUNCT__ "TSGetRHSJacobian"
1711 /*@C
1712    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1713 
1714    Not Collective, but parallel objects are returned if TS is parallel
1715 
1716    Input Parameter:
1717 .  ts  - The TS context obtained from TSCreate()
1718 
1719    Output Parameters:
1720 +  J   - The Jacobian J of F, where U_t = F(U,t)
1721 .  M   - The preconditioner matrix, usually the same as J
1722 - ctx - User-defined context for Jacobian evaluation routine
1723 
1724    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1725 
1726    Level: intermediate
1727 
1728 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
1729 
1730 .keywords: TS, timestep, get, matrix, Jacobian
1731 @*/
1732 PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1733 {
1734   PetscFunctionBegin;
1735   if (J) *J = ts->Arhs;
1736   if (M) *M = ts->B;
1737   if (ctx) *ctx = ts->jacP;
1738   PetscFunctionReturn(0);
1739 }
1740 
1741 #undef __FUNCT__
1742 #define __FUNCT__ "TSMonitorSolution"
1743 /*@C
1744    TSMonitorSolution - Monitors progress of the TS solvers by calling
1745    VecView() for the solution at each timestep
1746 
1747    Collective on TS
1748 
1749    Input Parameters:
1750 +  ts - the TS context
1751 .  step - current time-step
1752 .  ptime - current time
1753 -  dummy - either a viewer or PETSC_NULL
1754 
1755    Level: intermediate
1756 
1757 .keywords: TS,  vector, monitor, view
1758 
1759 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
1760 @*/
1761 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
1762 {
1763   PetscErrorCode ierr;
1764   PetscViewer    viewer = (PetscViewer) dummy;
1765 
1766   PetscFunctionBegin;
1767   if (!dummy) {
1768     viewer = PETSC_VIEWER_DRAW_(ts->comm);
1769   }
1770   ierr = VecView(x,viewer);CHKERRQ(ierr);
1771   PetscFunctionReturn(0);
1772 }
1773 
1774 
1775 
1776