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