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