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