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