xref: /petsc/src/ts/interface/ts.c (revision 5fa1c062e537a7dc06133d3274d2803bc931c085)
1 #define PETSCTS_DLL
2 
3 #include "src/ts/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_xmonitor - 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   PetscViewer    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","TSDefaultMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
97     if (flg) {
98       ierr = PetscViewerASCIIOpen(ts->comm,monfilename,&monviewer);CHKERRQ(ierr);
99       ierr = TSSetMonitor(ts,TSDefaultMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
100     }
101     ierr = PetscOptionsName("-ts_xmonitor","Monitor timestep size graphically","TSLGMonitor",&opt);CHKERRQ(ierr);
102     if (opt) {
103       ierr = TSSetMonitor(ts,TSLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
104     }
105     ierr = PetscOptionsName("-ts_vecmonitor","Monitor solution graphically","TSVecViewMonitor",&opt);CHKERRQ(ierr);
106     if (opt) {
107       ierr = TSSetMonitor(ts,TSVecViewMonitor,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->A,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->A,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->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
284       PetscStackPop;
285     }
286     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
287   }
288 
289   /* apply user-provided boundary conditions (only needed if these are time dependent) */
290   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
291   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
292 
293   PetscFunctionReturn(0);
294 }
295 
296 #undef __FUNCT__
297 #define __FUNCT__ "TSSetRHSFunction"
298 /*@C
299     TSSetRHSFunction - Sets the routine for evaluating the function,
300     F(t,u), where U_t = F(t,u).
301 
302     Collective on TS
303 
304     Input Parameters:
305 +   ts - the TS context obtained from TSCreate()
306 .   f - routine for evaluating the right-hand-side function
307 -   ctx - [optional] user-defined context for private data for the
308           function evaluation routine (may be PETSC_NULL)
309 
310     Calling sequence of func:
311 $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
312 
313 +   t - current timestep
314 .   u - input vector
315 .   F - function vector
316 -   ctx - [optional] user-defined function context
317 
318     Important:
319     The user MUST call either this routine or TSSetRHSMatrix().
320 
321     Level: beginner
322 
323 .keywords: TS, timestep, set, right-hand-side, function
324 
325 .seealso: TSSetRHSMatrix()
326 @*/
327 PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
328 {
329   PetscFunctionBegin;
330 
331   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
332   if (ts->problem_type == TS_LINEAR) {
333     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
334   }
335   ts->ops->rhsfunction = f;
336   ts->funP             = ctx;
337   PetscFunctionReturn(0);
338 }
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "TSSetRHSMatrix"
342 /*@C
343    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
344    Also sets the location to store A.
345 
346    Collective on TS
347 
348    Input Parameters:
349 +  ts  - the TS context obtained from TSCreate()
350 .  A   - matrix
351 .  B   - preconditioner matrix (usually same as A)
352 .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
353          if A is not a function of t.
354 -  ctx - [optional] user-defined context for private data for the
355           matrix evaluation routine (may be PETSC_NULL)
356 
357    Calling sequence of func:
358 $     func (TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
359 
360 +  t - current timestep
361 .  A - matrix A, where U_t = A(t) U
362 .  B - preconditioner matrix, usually the same as A
363 .  flag - flag indicating information about the preconditioner matrix
364           structure (same as flag in KSPSetOperators())
365 -  ctx - [optional] user-defined context for matrix evaluation routine
366 
367    Notes:
368    See KSPSetOperators() for important information about setting the flag
369    output parameter in the routine func().  Be sure to read this information!
370 
371    The routine func() takes Mat * as the matrix arguments rather than Mat.
372    This allows the matrix evaluation routine to replace A and/or B with a
373    completely new new matrix structure (not just different matrix elements)
374    when appropriate, for instance, if the nonzero structure is changing
375    throughout the global iterations.
376 
377    Important:
378    The user MUST call either this routine or TSSetRHSFunction().
379 
380    Level: beginner
381 
382 .keywords: TS, timestep, set, right-hand-side, matrix
383 
384 .seealso: TSSetRHSFunction()
385 @*/
386 PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSMatrix(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
387 {
388   PetscFunctionBegin;
389   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
390   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
391   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
392   PetscCheckSameComm(ts,1,A,2);
393   PetscCheckSameComm(ts,1,B,3);
394   if (ts->problem_type == TS_NONLINEAR) {
395     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
396   }
397 
398   ts->ops->rhsmatrix = f;
399   ts->jacP           = ctx;
400   ts->A              = A;
401   ts->B              = B;
402 
403   PetscFunctionReturn(0);
404 }
405 
406 #undef __FUNCT__
407 #define __FUNCT__ "TSSetRHSJacobian"
408 /*@C
409    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
410    where U_t = F(U,t), as well as the location to store the matrix.
411    Use TSSetRHSMatrix() for linear problems.
412 
413    Collective on TS
414 
415    Input Parameters:
416 +  ts  - the TS context obtained from TSCreate()
417 .  A   - Jacobian matrix
418 .  B   - preconditioner matrix (usually same as A)
419 .  f   - the Jacobian evaluation routine
420 -  ctx - [optional] user-defined context for private data for the
421          Jacobian evaluation routine (may be PETSC_NULL)
422 
423    Calling sequence of func:
424 $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
425 
426 +  t - current timestep
427 .  u - input vector
428 .  A - matrix A, where U_t = A(t)u
429 .  B - preconditioner matrix, usually the same as A
430 .  flag - flag indicating information about the preconditioner matrix
431           structure (same as flag in KSPSetOperators())
432 -  ctx - [optional] user-defined context for matrix evaluation routine
433 
434    Notes:
435    See KSPSetOperators() for important information about setting the flag
436    output parameter in the routine func().  Be sure to read this information!
437 
438    The routine func() takes Mat * as the matrix arguments rather than Mat.
439    This allows the matrix evaluation routine to replace A and/or B with a
440    completely new new matrix structure (not just different matrix elements)
441    when appropriate, for instance, if the nonzero structure is changing
442    throughout the global iterations.
443 
444    Level: beginner
445 
446 .keywords: TS, timestep, set, right-hand-side, Jacobian
447 
448 .seealso: TSDefaultComputeJacobianColor(),
449           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetRHSMatrix()
450 
451 @*/
452 PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
453 {
454   PetscFunctionBegin;
455   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
456   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
457   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
458   PetscCheckSameComm(ts,1,A,2);
459   PetscCheckSameComm(ts,1,B,3);
460   if (ts->problem_type != TS_NONLINEAR) {
461     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
462   }
463 
464   ts->ops->rhsjacobian = f;
465   ts->jacP             = ctx;
466   ts->A                = A;
467   ts->B                = B;
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "TSComputeRHSBoundaryConditions"
473 /*
474    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
475 
476    Note: If the user did not provide a function but merely a matrix,
477    this routine applies the matrix.
478 */
479 PetscErrorCode TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
480 {
481   PetscErrorCode ierr;
482 
483   PetscFunctionBegin;
484   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
485   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
486   PetscCheckSameComm(ts,1,x,3);
487 
488   if (ts->ops->rhsbc) {
489     PetscStackPush("TS user boundary condition function");
490     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
491     PetscStackPop;
492     PetscFunctionReturn(0);
493   }
494 
495   PetscFunctionReturn(0);
496 }
497 
498 #undef __FUNCT__
499 #define __FUNCT__ "TSSetRHSBoundaryConditions"
500 /*@C
501     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
502     boundary conditions for the function F.
503 
504     Collective on TS
505 
506     Input Parameters:
507 +   ts - the TS context obtained from TSCreate()
508 .   f - routine for evaluating the boundary condition function
509 -   ctx - [optional] user-defined context for private data for the
510           function evaluation routine (may be PETSC_NULL)
511 
512     Calling sequence of func:
513 $     func (TS ts,PetscReal t,Vec F,void *ctx);
514 
515 +   t - current timestep
516 .   F - function vector
517 -   ctx - [optional] user-defined function context
518 
519     Level: intermediate
520 
521 .keywords: TS, timestep, set, boundary conditions, function
522 @*/
523 PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSBoundaryConditions(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
524 {
525   PetscFunctionBegin;
526 
527   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
528   if (ts->problem_type != TS_LINEAR) {
529     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
530   }
531   ts->ops->rhsbc = f;
532   ts->bcP        = ctx;
533   PetscFunctionReturn(0);
534 }
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "TSView"
538 /*@C
539     TSView - Prints the TS data structure.
540 
541     Collective on TS
542 
543     Input Parameters:
544 +   ts - the TS context obtained from TSCreate()
545 -   viewer - visualization context
546 
547     Options Database Key:
548 .   -ts_view - calls TSView() at end of TSStep()
549 
550     Notes:
551     The available visualization contexts include
552 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
553 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
554          output where only the first processor opens
555          the file.  All other processors send their
556          data to the first processor to print.
557 
558     The user can open an alternative visualization context with
559     PetscViewerASCIIOpen() - output to a specified file.
560 
561     Level: beginner
562 
563 .keywords: TS, timestep, view
564 
565 .seealso: PetscViewerASCIIOpen()
566 @*/
567 PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer)
568 {
569   PetscErrorCode ierr;
570   char           *type;
571   PetscTruth     iascii,isstring;
572 
573   PetscFunctionBegin;
574   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
575   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
576   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
577   PetscCheckSameComm(ts,1,viewer,2);
578 
579   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
580   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
581   if (iascii) {
582     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
583     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
584     if (type) {
585       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
586     } else {
587       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
588     }
589     if (ts->ops->view) {
590       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
591       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
592       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
593     }
594     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
595     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
596     if (ts->problem_type == TS_NONLINEAR) {
597       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
598     }
599     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
600   } else if (isstring) {
601     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
602     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
603   }
604   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
605   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
606   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
607   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
608   PetscFunctionReturn(0);
609 }
610 
611 
612 #undef __FUNCT__
613 #define __FUNCT__ "TSSetApplicationContext"
614 /*@C
615    TSSetApplicationContext - Sets an optional user-defined context for
616    the timesteppers.
617 
618    Collective on TS
619 
620    Input Parameters:
621 +  ts - the TS context obtained from TSCreate()
622 -  usrP - optional user context
623 
624    Level: intermediate
625 
626 .keywords: TS, timestep, set, application, context
627 
628 .seealso: TSGetApplicationContext()
629 @*/
630 PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP)
631 {
632   PetscFunctionBegin;
633   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
634   ts->user = usrP;
635   PetscFunctionReturn(0);
636 }
637 
638 #undef __FUNCT__
639 #define __FUNCT__ "TSGetApplicationContext"
640 /*@C
641     TSGetApplicationContext - Gets the user-defined context for the
642     timestepper.
643 
644     Not Collective
645 
646     Input Parameter:
647 .   ts - the TS context obtained from TSCreate()
648 
649     Output Parameter:
650 .   usrP - user context
651 
652     Level: intermediate
653 
654 .keywords: TS, timestep, get, application, context
655 
656 .seealso: TSSetApplicationContext()
657 @*/
658 PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP)
659 {
660   PetscFunctionBegin;
661   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
662   *usrP = ts->user;
663   PetscFunctionReturn(0);
664 }
665 
666 #undef __FUNCT__
667 #define __FUNCT__ "TSGetTimeStepNumber"
668 /*@
669    TSGetTimeStepNumber - Gets the current number of timesteps.
670 
671    Not Collective
672 
673    Input Parameter:
674 .  ts - the TS context obtained from TSCreate()
675 
676    Output Parameter:
677 .  iter - number steps so far
678 
679    Level: intermediate
680 
681 .keywords: TS, timestep, get, iteration, number
682 @*/
683 PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter)
684 {
685   PetscFunctionBegin;
686   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
687   PetscValidIntPointer(iter,2);
688   *iter = ts->steps;
689   PetscFunctionReturn(0);
690 }
691 
692 #undef __FUNCT__
693 #define __FUNCT__ "TSSetInitialTimeStep"
694 /*@
695    TSSetInitialTimeStep - Sets the initial timestep to be used,
696    as well as the initial time.
697 
698    Collective on TS
699 
700    Input Parameters:
701 +  ts - the TS context obtained from TSCreate()
702 .  initial_time - the initial time
703 -  time_step - the size of the timestep
704 
705    Level: intermediate
706 
707 .seealso: TSSetTimeStep(), TSGetTimeStep()
708 
709 .keywords: TS, set, initial, timestep
710 @*/
711 PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
712 {
713   PetscFunctionBegin;
714   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
715   ts->time_step         = time_step;
716   ts->initial_time_step = time_step;
717   ts->ptime             = initial_time;
718   PetscFunctionReturn(0);
719 }
720 
721 #undef __FUNCT__
722 #define __FUNCT__ "TSSetTimeStep"
723 /*@
724    TSSetTimeStep - Allows one to reset the timestep at any time,
725    useful for simple pseudo-timestepping codes.
726 
727    Collective on TS
728 
729    Input Parameters:
730 +  ts - the TS context obtained from TSCreate()
731 -  time_step - the size of the timestep
732 
733    Level: intermediate
734 
735 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
736 
737 .keywords: TS, set, timestep
738 @*/
739 PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step)
740 {
741   PetscFunctionBegin;
742   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
743   ts->time_step = time_step;
744   PetscFunctionReturn(0);
745 }
746 
747 #undef __FUNCT__
748 #define __FUNCT__ "TSGetTimeStep"
749 /*@
750    TSGetTimeStep - Gets the current timestep size.
751 
752    Not Collective
753 
754    Input Parameter:
755 .  ts - the TS context obtained from TSCreate()
756 
757    Output Parameter:
758 .  dt - the current timestep size
759 
760    Level: intermediate
761 
762 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
763 
764 .keywords: TS, get, timestep
765 @*/
766 PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt)
767 {
768   PetscFunctionBegin;
769   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
770   PetscValidDoublePointer(dt,2);
771   *dt = ts->time_step;
772   PetscFunctionReturn(0);
773 }
774 
775 #undef __FUNCT__
776 #define __FUNCT__ "TSGetSolution"
777 /*@
778    TSGetSolution - Returns the solution at the present timestep. It
779    is valid to call this routine inside the function that you are evaluating
780    in order to move to the new timestep. This vector not changed until
781    the solution at the next timestep has been calculated.
782 
783    Not Collective, but Vec returned is parallel if TS is parallel
784 
785    Input Parameter:
786 .  ts - the TS context obtained from TSCreate()
787 
788    Output Parameter:
789 .  v - the vector containing the solution
790 
791    Level: intermediate
792 
793 .seealso: TSGetTimeStep()
794 
795 .keywords: TS, timestep, get, solution
796 @*/
797 PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v)
798 {
799   PetscFunctionBegin;
800   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
801   PetscValidPointer(v,2);
802   *v = ts->vec_sol_always;
803   PetscFunctionReturn(0);
804 }
805 
806 /* ----- Routines to initialize and destroy a timestepper ---- */
807 #undef __FUNCT__
808 #define __FUNCT__ "TSSetProblemType"
809 /*@
810   TSSetProblemType - Sets the type of problem to be solved.
811 
812   Not collective
813 
814   Input Parameters:
815 + ts   - The TS
816 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
817 .vb
818          U_t = A U
819          U_t = A(t) U
820          U_t = F(t,U)
821 .ve
822 
823    Level: beginner
824 
825 .keywords: TS, problem type
826 .seealso: TSSetUp(), TSProblemType, TS
827 @*/
828 PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type)
829 {
830   PetscFunctionBegin;
831   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
832   ts->problem_type = type;
833   PetscFunctionReturn(0);
834 }
835 
836 #undef __FUNCT__
837 #define __FUNCT__ "TSGetProblemType"
838 /*@C
839   TSGetProblemType - Gets the type of problem to be solved.
840 
841   Not collective
842 
843   Input Parameter:
844 . ts   - The TS
845 
846   Output Parameter:
847 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
848 .vb
849          U_t = A U
850          U_t = A(t) U
851          U_t = F(t,U)
852 .ve
853 
854    Level: beginner
855 
856 .keywords: TS, problem type
857 .seealso: TSSetUp(), TSProblemType, TS
858 @*/
859 PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type)
860 {
861   PetscFunctionBegin;
862   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
863   PetscValidIntPointer(type,2);
864   *type = ts->problem_type;
865   PetscFunctionReturn(0);
866 }
867 
868 #undef __FUNCT__
869 #define __FUNCT__ "TSSetUp"
870 /*@
871    TSSetUp - Sets up the internal data structures for the later use
872    of a timestepper.
873 
874    Collective on TS
875 
876    Input Parameter:
877 .  ts - the TS context obtained from TSCreate()
878 
879    Notes:
880    For basic use of the TS solvers the user need not explicitly call
881    TSSetUp(), since these actions will automatically occur during
882    the call to TSStep().  However, if one wishes to control this
883    phase separately, TSSetUp() should be called after TSCreate()
884    and optional routines of the form TSSetXXX(), but before TSStep().
885 
886    Level: advanced
887 
888 .keywords: TS, timestep, setup
889 
890 .seealso: TSCreate(), TSStep(), TSDestroy()
891 @*/
892 PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts)
893 {
894   PetscErrorCode ierr;
895 
896   PetscFunctionBegin;
897   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
898   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
899   if (!ts->type_name) {
900     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
901   }
902   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
903   ts->setupcalled = 1;
904   PetscFunctionReturn(0);
905 }
906 
907 #undef __FUNCT__
908 #define __FUNCT__ "TSDestroy"
909 /*@
910    TSDestroy - Destroys the timestepper context that was created
911    with TSCreate().
912 
913    Collective on TS
914 
915    Input Parameter:
916 .  ts - the TS context obtained from TSCreate()
917 
918    Level: beginner
919 
920 .keywords: TS, timestepper, destroy
921 
922 .seealso: TSCreate(), TSSetUp(), TSSolve()
923 @*/
924 PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts)
925 {
926   PetscErrorCode ierr;
927 
928   PetscFunctionBegin;
929   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
930   if (--ts->refct > 0) PetscFunctionReturn(0);
931 
932   /* if memory was published with AMS then destroy it */
933   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
934 
935   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
936   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
937   ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);
938   ierr = TSClearMonitor(ts);CHKERRQ(ierr);
939   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
940   PetscFunctionReturn(0);
941 }
942 
943 #undef __FUNCT__
944 #define __FUNCT__ "TSGetSNES"
945 /*@
946    TSGetSNES - Returns the SNES (nonlinear solver) associated with
947    a TS (timestepper) context. Valid only for nonlinear problems.
948 
949    Not Collective, but SNES is parallel if TS is parallel
950 
951    Input Parameter:
952 .  ts - the TS context obtained from TSCreate()
953 
954    Output Parameter:
955 .  snes - the nonlinear solver context
956 
957    Notes:
958    The user can then directly manipulate the SNES context to set various
959    options, etc.  Likewise, the user can then extract and manipulate the
960    KSP, KSP, and PC contexts as well.
961 
962    TSGetSNES() does not work for integrators that do not use SNES; in
963    this case TSGetSNES() returns PETSC_NULL in snes.
964 
965    Level: beginner
966 
967 .keywords: timestep, get, SNES
968 @*/
969 PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes)
970 {
971   PetscFunctionBegin;
972   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
973   PetscValidPointer(snes,2);
974   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
975   *snes = ts->snes;
976   PetscFunctionReturn(0);
977 }
978 
979 #undef __FUNCT__
980 #define __FUNCT__ "TSGetKSP"
981 /*@
982    TSGetKSP - Returns the KSP (linear solver) associated with
983    a TS (timestepper) context.
984 
985    Not Collective, but KSP is parallel if TS is parallel
986 
987    Input Parameter:
988 .  ts - the TS context obtained from TSCreate()
989 
990    Output Parameter:
991 .  ksp - the nonlinear solver context
992 
993    Notes:
994    The user can then directly manipulate the KSP context to set various
995    options, etc.  Likewise, the user can then extract and manipulate the
996    KSP and PC contexts as well.
997 
998    TSGetKSP() does not work for integrators that do not use KSP;
999    in this case TSGetKSP() returns PETSC_NULL in ksp.
1000 
1001    Level: beginner
1002 
1003 .keywords: timestep, get, KSP
1004 @*/
1005 PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp)
1006 {
1007   PetscFunctionBegin;
1008   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1009   PetscValidPointer(ksp,2);
1010   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1011   *ksp = ts->ksp;
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 /* ----------- Routines to set solver parameters ---------- */
1016 
1017 #undef __FUNCT__
1018 #define __FUNCT__ "TSGetDuration"
1019 /*@
1020    TSGetDuration - Gets the maximum number of timesteps to use and
1021    maximum time for iteration.
1022 
1023    Collective on TS
1024 
1025    Input Parameters:
1026 +  ts       - the TS context obtained from TSCreate()
1027 .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1028 -  maxtime  - final time to iterate to, or PETSC_NULL
1029 
1030    Level: intermediate
1031 
1032 .keywords: TS, timestep, get, maximum, iterations, time
1033 @*/
1034 PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1035 {
1036   PetscFunctionBegin;
1037   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1038   if (maxsteps) {
1039     PetscValidIntPointer(maxsteps,2);
1040     *maxsteps = ts->max_steps;
1041   }
1042   if (maxtime ) {
1043     PetscValidScalarPointer(maxtime,3);
1044     *maxtime  = ts->max_time;
1045   }
1046   PetscFunctionReturn(0);
1047 }
1048 
1049 #undef __FUNCT__
1050 #define __FUNCT__ "TSSetDuration"
1051 /*@
1052    TSSetDuration - Sets the maximum number of timesteps to use and
1053    maximum time for iteration.
1054 
1055    Collective on TS
1056 
1057    Input Parameters:
1058 +  ts - the TS context obtained from TSCreate()
1059 .  maxsteps - maximum number of iterations to use
1060 -  maxtime - final time to iterate to
1061 
1062    Options Database Keys:
1063 .  -ts_max_steps <maxsteps> - Sets maxsteps
1064 .  -ts_max_time <maxtime> - Sets maxtime
1065 
1066    Notes:
1067    The default maximum number of iterations is 5000. Default time is 5.0
1068 
1069    Level: intermediate
1070 
1071 .keywords: TS, timestep, set, maximum, iterations
1072 @*/
1073 PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1074 {
1075   PetscFunctionBegin;
1076   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1077   ts->max_steps = maxsteps;
1078   ts->max_time  = maxtime;
1079   PetscFunctionReturn(0);
1080 }
1081 
1082 #undef __FUNCT__
1083 #define __FUNCT__ "TSSetSolution"
1084 /*@
1085    TSSetSolution - Sets the initial solution vector
1086    for use by the TS routines.
1087 
1088    Collective on TS and Vec
1089 
1090    Input Parameters:
1091 +  ts - the TS context obtained from TSCreate()
1092 -  x - the solution vector
1093 
1094    Level: beginner
1095 
1096 .keywords: TS, timestep, set, solution, initial conditions
1097 @*/
1098 PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x)
1099 {
1100   PetscFunctionBegin;
1101   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1102   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1103   ts->vec_sol        = ts->vec_sol_always = x;
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 #undef __FUNCT__
1108 #define __FUNCT__ "TSDefaultRhsBC"
1109 /*@
1110   TSDefaultRhsBC - The default boundary condition function which does nothing.
1111 
1112   Collective on TS
1113 
1114   Input Parameters:
1115 + ts  - The TS context obtained from TSCreate()
1116 . rhs - The Rhs
1117 - ctx - The user-context
1118 
1119   Level: developer
1120 
1121 .keywords: TS, Rhs, boundary conditions
1122 @*/
1123 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1124 {
1125   PetscFunctionBegin;
1126   PetscFunctionReturn(0);
1127 }
1128 
1129 #undef __FUNCT__
1130 #define __FUNCT__ "TSSetSystemMatrixBC"
1131 /*@C
1132   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1133   to the system matrix and preconditioner of each system.
1134 
1135   Collective on TS
1136 
1137   Input Parameters:
1138 + ts   - The TS context obtained from TSCreate()
1139 - func - The function
1140 
1141   Calling sequence of func:
1142 . func (TS ts, Mat A, Mat B, void *ctx);
1143 
1144 + A   - The current system matrix
1145 . B   - The current preconditioner
1146 - ctx - The user-context
1147 
1148   Level: intermediate
1149 
1150 .keywords: TS, System matrix, boundary conditions
1151 @*/
1152 PetscErrorCode PETSCTS_DLLEXPORT TSSetSystemMatrixBC(TS ts, PetscErrorCode (*func)(TS, Mat, Mat, void *))
1153 {
1154   PetscFunctionBegin;
1155   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1156   ts->ops->applymatrixbc = func;
1157   PetscFunctionReturn(0);
1158 }
1159 
1160 #undef __FUNCT__
1161 #define __FUNCT__ "TSDefaultSystemMatrixBC"
1162 /*@
1163   TSDefaultSystemMatrixBC - The default boundary condition function which
1164   does nothing.
1165 
1166   Collective on TS
1167 
1168   Input Parameters:
1169 + ts  - The TS context obtained from TSCreate()
1170 . A   - The system matrix
1171 . B   - The preconditioner
1172 - ctx - The user-context
1173 
1174   Level: developer
1175 
1176 .keywords: TS, System matrix, boundary conditions
1177 @*/
1178 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1179 {
1180   PetscFunctionBegin;
1181   PetscFunctionReturn(0);
1182 }
1183 
1184 #undef __FUNCT__
1185 #define __FUNCT__ "TSSetSolutionBC"
1186 /*@C
1187   TSSetSolutionBC - Sets the function which applies boundary conditions
1188   to the solution of each system. This is necessary in nonlinear systems
1189   which time dependent boundary conditions.
1190 
1191   Collective on TS
1192 
1193   Input Parameters:
1194 + ts   - The TS context obtained from TSCreate()
1195 - func - The function
1196 
1197   Calling sequence of func:
1198 . func (TS ts, Vec rsol, void *ctx);
1199 
1200 + sol - The current solution vector
1201 - ctx - The user-context
1202 
1203   Level: intermediate
1204 
1205 .keywords: TS, solution, boundary conditions
1206 @*/
1207 PetscErrorCode PETSCTS_DLLEXPORT TSSetSolutionBC(TS ts, PetscErrorCode (*func)(TS, Vec, void *))
1208 {
1209   PetscFunctionBegin;
1210   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1211   ts->ops->applysolbc = func;
1212   PetscFunctionReturn(0);
1213 }
1214 
1215 #undef __FUNCT__
1216 #define __FUNCT__ "TSDefaultSolutionBC"
1217 /*@
1218   TSDefaultSolutionBC - The default boundary condition function which
1219   does nothing.
1220 
1221   Collective on TS
1222 
1223   Input Parameters:
1224 + ts  - The TS context obtained from TSCreate()
1225 . sol - The solution
1226 - ctx - The user-context
1227 
1228   Level: developer
1229 
1230 .keywords: TS, solution, boundary conditions
1231 @*/
1232 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1233 {
1234   PetscFunctionBegin;
1235   PetscFunctionReturn(0);
1236 }
1237 
1238 #undef __FUNCT__
1239 #define __FUNCT__ "TSSetPreStep"
1240 /*@C
1241   TSSetPreStep - Sets the general-purpose function
1242   called once at the beginning of time stepping.
1243 
1244   Collective on TS
1245 
1246   Input Parameters:
1247 + ts   - The TS context obtained from TSCreate()
1248 - func - The function
1249 
1250   Calling sequence of func:
1251 . func (TS ts);
1252 
1253   Level: intermediate
1254 
1255 .keywords: TS, timestep
1256 @*/
1257 PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1258 {
1259   PetscFunctionBegin;
1260   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1261   ts->ops->prestep = func;
1262   PetscFunctionReturn(0);
1263 }
1264 
1265 #undef __FUNCT__
1266 #define __FUNCT__ "TSDefaultPreStep"
1267 /*@
1268   TSDefaultPreStep - The default pre-stepping function which does nothing.
1269 
1270   Collective on TS
1271 
1272   Input Parameters:
1273 . ts  - The TS context obtained from TSCreate()
1274 
1275   Level: developer
1276 
1277 .keywords: TS, timestep
1278 @*/
1279 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts)
1280 {
1281   PetscFunctionBegin;
1282   PetscFunctionReturn(0);
1283 }
1284 
1285 #undef __FUNCT__
1286 #define __FUNCT__ "TSSetUpdate"
1287 /*@C
1288   TSSetUpdate - Sets the general-purpose update function called
1289   at the beginning of every time step. This function can change
1290   the time step.
1291 
1292   Collective on TS
1293 
1294   Input Parameters:
1295 + ts   - The TS context obtained from TSCreate()
1296 - func - The function
1297 
1298   Calling sequence of func:
1299 . func (TS ts, double t, double *dt);
1300 
1301 + t   - The current time
1302 - dt  - The current time step
1303 
1304   Level: intermediate
1305 
1306 .keywords: TS, update, timestep
1307 @*/
1308 PetscErrorCode PETSCTS_DLLEXPORT TSSetUpdate(TS ts, PetscErrorCode (*func)(TS, PetscReal, PetscReal *))
1309 {
1310   PetscFunctionBegin;
1311   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1312   ts->ops->update = func;
1313   PetscFunctionReturn(0);
1314 }
1315 
1316 #undef __FUNCT__
1317 #define __FUNCT__ "TSDefaultUpdate"
1318 /*@
1319   TSDefaultUpdate - The default update function which does nothing.
1320 
1321   Collective on TS
1322 
1323   Input Parameters:
1324 + ts  - The TS context obtained from TSCreate()
1325 - t   - The current time
1326 
1327   Output Parameters:
1328 . dt  - The current time step
1329 
1330   Level: developer
1331 
1332 .keywords: TS, update, timestep
1333 @*/
1334 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1335 {
1336   PetscFunctionBegin;
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 #undef __FUNCT__
1341 #define __FUNCT__ "TSSetPostStep"
1342 /*@C
1343   TSSetPostStep - Sets the general-purpose function
1344   called once at the end of time stepping.
1345 
1346   Collective on TS
1347 
1348   Input Parameters:
1349 + ts   - The TS context obtained from TSCreate()
1350 - func - The function
1351 
1352   Calling sequence of func:
1353 . func (TS ts);
1354 
1355   Level: intermediate
1356 
1357 .keywords: TS, timestep
1358 @*/
1359 PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1360 {
1361   PetscFunctionBegin;
1362   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1363   ts->ops->poststep = func;
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 #undef __FUNCT__
1368 #define __FUNCT__ "TSDefaultPostStep"
1369 /*@
1370   TSDefaultPostStep - The default post-stepping function which does nothing.
1371 
1372   Collective on TS
1373 
1374   Input Parameters:
1375 . ts  - The TS context obtained from TSCreate()
1376 
1377   Level: developer
1378 
1379 .keywords: TS, timestep
1380 @*/
1381 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts)
1382 {
1383   PetscFunctionBegin;
1384   PetscFunctionReturn(0);
1385 }
1386 
1387 /* ------------ Routines to set performance monitoring options ----------- */
1388 
1389 #undef __FUNCT__
1390 #define __FUNCT__ "TSSetMonitor"
1391 /*@C
1392    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1393    timestep to display the iteration's  progress.
1394 
1395    Collective on TS
1396 
1397    Input Parameters:
1398 +  ts - the TS context obtained from TSCreate()
1399 .  func - monitoring routine
1400 .  mctx - [optional] user-defined context for private data for the
1401              monitor routine (use PETSC_NULL if no context is desired)
1402 -  monitordestroy - [optional] routine that frees monitor context
1403           (may be PETSC_NULL)
1404 
1405    Calling sequence of func:
1406 $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1407 
1408 +    ts - the TS context
1409 .    steps - iteration number
1410 .    time - current time
1411 .    x - current iterate
1412 -    mctx - [optional] monitoring context
1413 
1414    Notes:
1415    This routine adds an additional monitor to the list of monitors that
1416    already has been loaded.
1417 
1418    Level: intermediate
1419 
1420 .keywords: TS, timestep, set, monitor
1421 
1422 .seealso: TSDefaultMonitor(), TSClearMonitor()
1423 @*/
1424 PetscErrorCode PETSCTS_DLLEXPORT TSSetMonitor(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1425 {
1426   PetscFunctionBegin;
1427   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1428   if (ts->numbermonitors >= MAXTSMONITORS) {
1429     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1430   }
1431   ts->monitor[ts->numbermonitors]           = monitor;
1432   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1433   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1434   PetscFunctionReturn(0);
1435 }
1436 
1437 #undef __FUNCT__
1438 #define __FUNCT__ "TSClearMonitor"
1439 /*@C
1440    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1441 
1442    Collective on TS
1443 
1444    Input Parameters:
1445 .  ts - the TS context obtained from TSCreate()
1446 
1447    Notes:
1448    There is no way to remove a single, specific monitor.
1449 
1450    Level: intermediate
1451 
1452 .keywords: TS, timestep, set, monitor
1453 
1454 .seealso: TSDefaultMonitor(), TSSetMonitor()
1455 @*/
1456 PetscErrorCode PETSCTS_DLLEXPORT TSClearMonitor(TS ts)
1457 {
1458   PetscErrorCode ierr;
1459   PetscInt       i;
1460 
1461   PetscFunctionBegin;
1462   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1463   for (i=0; i<ts->numbermonitors; i++) {
1464     if (ts->mdestroy[i]) {
1465       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
1466     }
1467   }
1468   ts->numbermonitors = 0;
1469   PetscFunctionReturn(0);
1470 }
1471 
1472 #undef __FUNCT__
1473 #define __FUNCT__ "TSDefaultMonitor"
1474 /*@
1475    TSDefaultMonitor - Sets the Default monitor
1476 @*/
1477 PetscErrorCode TSDefaultMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1478 {
1479   PetscErrorCode ierr;
1480   PetscViewer    viewer = (PetscViewer)ctx;
1481 
1482   PetscFunctionBegin;
1483   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
1484   ierr = PetscViewerASCIIPrintf(viewer,"timestep %D dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1485   PetscFunctionReturn(0);
1486 }
1487 
1488 #undef __FUNCT__
1489 #define __FUNCT__ "TSStep"
1490 /*@
1491    TSStep - Steps the requested number of timesteps.
1492 
1493    Collective on TS
1494 
1495    Input Parameter:
1496 .  ts - the TS context obtained from TSCreate()
1497 
1498    Output Parameters:
1499 +  steps - number of iterations until termination
1500 -  ptime - time until termination
1501 
1502    Level: beginner
1503 
1504 .keywords: TS, timestep, solve
1505 
1506 .seealso: TSCreate(), TSSetUp(), TSDestroy()
1507 @*/
1508 PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1509 {
1510   PetscErrorCode ierr;
1511 
1512   PetscFunctionBegin;
1513   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1514   if (!ts->setupcalled) {
1515     ierr = TSSetUp(ts);CHKERRQ(ierr);
1516   }
1517 
1518   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1519   ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1520   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1521   ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1522   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1523 
1524   if (!PetscPreLoadingOn) {
1525     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1526   }
1527   PetscFunctionReturn(0);
1528 }
1529 
1530 #undef __FUNCT__
1531 #define __FUNCT__ "TSMonitor"
1532 /*
1533      Runs the user provided monitor routines, if they exists.
1534 */
1535 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1536 {
1537   PetscErrorCode ierr;
1538   PetscInt i,n = ts->numbermonitors;
1539 
1540   PetscFunctionBegin;
1541   for (i=0; i<n; i++) {
1542     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1543   }
1544   PetscFunctionReturn(0);
1545 }
1546 
1547 /* ------------------------------------------------------------------------*/
1548 
1549 #undef __FUNCT__
1550 #define __FUNCT__ "TSLGMonitorCreate"
1551 /*@C
1552    TSLGMonitorCreate - Creates a line graph context for use with
1553    TS to monitor convergence of preconditioned residual norms.
1554 
1555    Collective on TS
1556 
1557    Input Parameters:
1558 +  host - the X display to open, or null for the local machine
1559 .  label - the title to put in the title bar
1560 .  x, y - the screen coordinates of the upper left coordinate of the window
1561 -  m, n - the screen width and height in pixels
1562 
1563    Output Parameter:
1564 .  draw - the drawing context
1565 
1566    Options Database Key:
1567 .  -ts_xmonitor - automatically sets line graph monitor
1568 
1569    Notes:
1570    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1571 
1572    Level: intermediate
1573 
1574 .keywords: TS, monitor, line graph, residual, seealso
1575 
1576 .seealso: TSLGMonitorDestroy(), TSSetMonitor()
1577 
1578 @*/
1579 PetscErrorCode PETSCTS_DLLEXPORT TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1580 {
1581   PetscDraw      win;
1582   PetscErrorCode ierr;
1583 
1584   PetscFunctionBegin;
1585   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1586   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1587   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1588   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1589 
1590   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1591   PetscFunctionReturn(0);
1592 }
1593 
1594 #undef __FUNCT__
1595 #define __FUNCT__ "TSLGMonitor"
1596 PetscErrorCode TSLGMonitor(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1597 {
1598   PetscDrawLG    lg = (PetscDrawLG) monctx;
1599   PetscReal      x,y = ptime;
1600   PetscErrorCode ierr;
1601 
1602   PetscFunctionBegin;
1603   if (!monctx) {
1604     MPI_Comm    comm;
1605     PetscViewer viewer;
1606 
1607     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1608     viewer = PETSC_VIEWER_DRAW_(comm);
1609     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
1610   }
1611 
1612   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1613   x = (PetscReal)n;
1614   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1615   if (n < 20 || (n % 5)) {
1616     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1617   }
1618   PetscFunctionReturn(0);
1619 }
1620 
1621 #undef __FUNCT__
1622 #define __FUNCT__ "TSLGMonitorDestroy"
1623 /*@C
1624    TSLGMonitorDestroy - Destroys a line graph context that was created
1625    with TSLGMonitorCreate().
1626 
1627    Collective on PetscDrawLG
1628 
1629    Input Parameter:
1630 .  draw - the drawing context
1631 
1632    Level: intermediate
1633 
1634 .keywords: TS, monitor, line graph, destroy
1635 
1636 .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1637 @*/
1638 PetscErrorCode PETSCTS_DLLEXPORT TSLGMonitorDestroy(PetscDrawLG drawlg)
1639 {
1640   PetscDraw      draw;
1641   PetscErrorCode ierr;
1642 
1643   PetscFunctionBegin;
1644   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1645   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1646   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1647   PetscFunctionReturn(0);
1648 }
1649 
1650 #undef __FUNCT__
1651 #define __FUNCT__ "TSGetTime"
1652 /*@
1653    TSGetTime - Gets the current time.
1654 
1655    Not Collective
1656 
1657    Input Parameter:
1658 .  ts - the TS context obtained from TSCreate()
1659 
1660    Output Parameter:
1661 .  t  - the current time
1662 
1663    Contributed by: Matthew Knepley
1664 
1665    Level: beginner
1666 
1667 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1668 
1669 .keywords: TS, get, time
1670 @*/
1671 PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t)
1672 {
1673   PetscFunctionBegin;
1674   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1675   PetscValidDoublePointer(t,2);
1676   *t = ts->ptime;
1677   PetscFunctionReturn(0);
1678 }
1679 
1680 #undef __FUNCT__
1681 #define __FUNCT__ "TSSetOptionsPrefix"
1682 /*@C
1683    TSSetOptionsPrefix - Sets the prefix used for searching for all
1684    TS options in the database.
1685 
1686    Collective on TS
1687 
1688    Input Parameter:
1689 +  ts     - The TS context
1690 -  prefix - The prefix to prepend to all option names
1691 
1692    Notes:
1693    A hyphen (-) must NOT be given at the beginning of the prefix name.
1694    The first character of all runtime options is AUTOMATICALLY the
1695    hyphen.
1696 
1697    Contributed by: Matthew Knepley
1698 
1699    Level: advanced
1700 
1701 .keywords: TS, set, options, prefix, database
1702 
1703 .seealso: TSSetFromOptions()
1704 
1705 @*/
1706 PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[])
1707 {
1708   PetscErrorCode ierr;
1709 
1710   PetscFunctionBegin;
1711   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1712   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1713   switch(ts->problem_type) {
1714     case TS_NONLINEAR:
1715       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1716       break;
1717     case TS_LINEAR:
1718       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1719       break;
1720   }
1721   PetscFunctionReturn(0);
1722 }
1723 
1724 
1725 #undef __FUNCT__
1726 #define __FUNCT__ "TSAppendOptionsPrefix"
1727 /*@C
1728    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1729    TS options in the database.
1730 
1731    Collective on TS
1732 
1733    Input Parameter:
1734 +  ts     - The TS context
1735 -  prefix - The prefix to prepend to all option names
1736 
1737    Notes:
1738    A hyphen (-) must NOT be given at the beginning of the prefix name.
1739    The first character of all runtime options is AUTOMATICALLY the
1740    hyphen.
1741 
1742    Contributed by: Matthew Knepley
1743 
1744    Level: advanced
1745 
1746 .keywords: TS, append, options, prefix, database
1747 
1748 .seealso: TSGetOptionsPrefix()
1749 
1750 @*/
1751 PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[])
1752 {
1753   PetscErrorCode ierr;
1754 
1755   PetscFunctionBegin;
1756   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1757   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1758   switch(ts->problem_type) {
1759     case TS_NONLINEAR:
1760       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1761       break;
1762     case TS_LINEAR:
1763       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1764       break;
1765   }
1766   PetscFunctionReturn(0);
1767 }
1768 
1769 #undef __FUNCT__
1770 #define __FUNCT__ "TSGetOptionsPrefix"
1771 /*@C
1772    TSGetOptionsPrefix - Sets the prefix used for searching for all
1773    TS options in the database.
1774 
1775    Not Collective
1776 
1777    Input Parameter:
1778 .  ts - The TS context
1779 
1780    Output Parameter:
1781 .  prefix - A pointer to the prefix string used
1782 
1783    Contributed by: Matthew Knepley
1784 
1785    Notes: On the fortran side, the user should pass in a string 'prifix' of
1786    sufficient length to hold the prefix.
1787 
1788    Level: intermediate
1789 
1790 .keywords: TS, get, options, prefix, database
1791 
1792 .seealso: TSAppendOptionsPrefix()
1793 @*/
1794 PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[])
1795 {
1796   PetscErrorCode ierr;
1797 
1798   PetscFunctionBegin;
1799   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1800   PetscValidPointer(prefix,2);
1801   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1802   PetscFunctionReturn(0);
1803 }
1804 
1805 #undef __FUNCT__
1806 #define __FUNCT__ "TSGetRHSMatrix"
1807 /*@C
1808    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1809 
1810    Not Collective, but parallel objects are returned if TS is parallel
1811 
1812    Input Parameter:
1813 .  ts  - The TS context obtained from TSCreate()
1814 
1815    Output Parameters:
1816 +  A   - The matrix A, where U_t = A(t) U
1817 .  M   - The preconditioner matrix, usually the same as A
1818 -  ctx - User-defined context for matrix evaluation routine
1819 
1820    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1821 
1822    Contributed by: Matthew Knepley
1823 
1824    Level: intermediate
1825 
1826 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1827 
1828 .keywords: TS, timestep, get, matrix
1829 
1830 @*/
1831 PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1832 {
1833   PetscFunctionBegin;
1834   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1835   if (A)   *A = ts->A;
1836   if (M)   *M = ts->B;
1837   if (ctx) *ctx = ts->jacP;
1838   PetscFunctionReturn(0);
1839 }
1840 
1841 #undef __FUNCT__
1842 #define __FUNCT__ "TSGetRHSJacobian"
1843 /*@C
1844    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1845 
1846    Not Collective, but parallel objects are returned if TS is parallel
1847 
1848    Input Parameter:
1849 .  ts  - The TS context obtained from TSCreate()
1850 
1851    Output Parameters:
1852 +  J   - The Jacobian J of F, where U_t = F(U,t)
1853 .  M   - The preconditioner matrix, usually the same as J
1854 - ctx - User-defined context for Jacobian evaluation routine
1855 
1856    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1857 
1858    Contributed by: Matthew Knepley
1859 
1860    Level: intermediate
1861 
1862 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1863 
1864 .keywords: TS, timestep, get, matrix, Jacobian
1865 @*/
1866 PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1867 {
1868   PetscErrorCode ierr;
1869 
1870   PetscFunctionBegin;
1871   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1872   PetscFunctionReturn(0);
1873 }
1874 
1875 #undef __FUNCT__
1876 #define __FUNCT__ "TSVecViewMonitor"
1877 /*@C
1878    TSVecViewMonitor - Monitors progress of the TS solvers by calling
1879    VecView() for the solution at each timestep
1880 
1881    Collective on TS
1882 
1883    Input Parameters:
1884 +  ts - the TS context
1885 .  step - current time-step
1886 .  ptime - current time
1887 -  dummy - either a viewer or PETSC_NULL
1888 
1889    Level: intermediate
1890 
1891 .keywords: TS,  vector, monitor, view
1892 
1893 .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
1894 @*/
1895 PetscErrorCode PETSCTS_DLLEXPORT TSVecViewMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
1896 {
1897   PetscErrorCode ierr;
1898   PetscViewer    viewer = (PetscViewer) dummy;
1899 
1900   PetscFunctionBegin;
1901   if (!viewer) {
1902     MPI_Comm comm;
1903     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1904     viewer = PETSC_VIEWER_DRAW_(comm);
1905   }
1906   ierr = VecView(x,viewer);CHKERRQ(ierr);
1907   PetscFunctionReturn(0);
1908 }
1909 
1910 
1911 
1912