xref: /petsc/src/ts/interface/ts.c (revision 6831982abb6453c2f3e25efecb5d0951d333371e)
1 #ifdef PETSC_RCS_HEADER
2 
3 #endif
4 
5 #include "src/ts/tsimpl.h"        /*I "ts.h"  I*/
6 
7 #undef __FUNC__
8 #define __FUNC__ "TSComputeRHSFunction"
9 /*
10    TSComputeRHSFunction - Evaluates the right-hand-side function.
11 
12    Note: If the user did not provide a function but merely a matrix,
13    this routine applies the matrix.
14 */
15 int TSComputeRHSFunction(TS ts,double t,Vec x, Vec y)
16 {
17   int ierr;
18 
19   PetscFunctionBegin;
20   PetscValidHeaderSpecific(ts,TS_COOKIE);
21   PetscValidHeader(x);  PetscValidHeader(y);
22 
23   if (ts->rhsfunction) {
24     PetscStackPush("TS user right-hand-side function");
25     ierr = (*ts->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
26     PetscStackPop;
27     PetscFunctionReturn(0);
28   }
29 
30   if (ts->rhsmatrix) { /* assemble matrix for this timestep */
31     MatStructure flg;
32     PetscStackPush("TS user right-hand-side matrix function");
33     ierr = (*ts->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
34     PetscStackPop;
35   }
36   ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
37 
38   /* apply user-provided boundary conditions (only needed if these are time dependent) */
39   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
40 
41   PetscFunctionReturn(0);
42 }
43 
44 #undef __FUNC__
45 #define __FUNC__ "TSSetRHSFunction"
46 /*@C
47     TSSetRHSFunction - Sets the routine for evaluating the function,
48     F(t,u), where U_t = F(t,u).
49 
50     Collective on TS
51 
52     Input Parameters:
53 +   ts - the TS context obtained from TSCreate()
54 .   f - routine for evaluating the right-hand-side function
55 -   ctx - [optional] user-defined context for private data for the
56           function evaluation routine (may be PETSC_NULL)
57 
58     Calling sequence of func:
59 $     func (TS ts,double t,Vec u,Vec F,void *ctx);
60 
61 +   t - current timestep
62 .   u - input vector
63 .   F - function vector
64 -   ctx - [optional] user-defined function context
65 
66     Important:
67     The user MUST call either this routine or TSSetRHSMatrix().
68 
69     Level: beginner
70 
71 .keywords: TS, timestep, set, right-hand-side, function
72 
73 .seealso: TSSetRHSMatrix()
74 @*/
75 int TSSetRHSFunction(TS ts,int (*f)(TS,double,Vec,Vec,void*),void *ctx)
76 {
77   PetscFunctionBegin;
78 
79   PetscValidHeaderSpecific(ts,TS_COOKIE);
80   if (ts->problem_type == TS_LINEAR) {
81     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Cannot set function for linear problem");
82   }
83   ts->rhsfunction = f;
84   ts->funP        = ctx;
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNC__
89 #define __FUNC__ "TSSetRHSMatrix"
90 /*@C
91    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
92    Also sets the location to store A.
93 
94    Collective on TS
95 
96    Input Parameters:
97 +  ts  - the TS context obtained from TSCreate()
98 .  A   - matrix
99 .  B   - preconditioner matrix (usually same as A)
100 .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
101          if A is not a function of t.
102 -  ctx - [optional] user-defined context for private data for the
103           matrix evaluation routine (may be PETSC_NULL)
104 
105    Calling sequence of func:
106 $     func (TS ts,double t,Mat *A,Mat *B,int *flag,void *ctx);
107 
108 +  t - current timestep
109 .  A - matrix A, where U_t = A(t) U
110 .  B - preconditioner matrix, usually the same as A
111 .  flag - flag indicating information about the preconditioner matrix
112           structure (same as flag in SLESSetOperators())
113 -  ctx - [optional] user-defined context for matrix evaluation routine
114 
115    Notes:
116    See SLESSetOperators() for important information about setting the flag
117    output parameter in the routine func().  Be sure to read this information!
118 
119    The routine func() takes Mat * as the matrix arguments rather than Mat.
120    This allows the matrix evaluation routine to replace A and/or B with a
121    completely new new matrix structure (not just different matrix elements)
122    when appropriate, for instance, if the nonzero structure is changing
123    throughout the global iterations.
124 
125    Important:
126    The user MUST call either this routine or TSSetRHSFunction().
127 
128    Level: beginner
129 
130 .keywords: TS, timestep, set, right-hand-side, matrix
131 
132 .seealso: TSSetRHSFunction()
133 @*/
134 int TSSetRHSMatrix(TS ts,Mat A, Mat B,int (*f)(TS,double,Mat*,Mat*,MatStructure*,void*),void *ctx)
135 {
136   PetscFunctionBegin;
137   PetscValidHeaderSpecific(ts,TS_COOKIE);
138   PetscValidHeaderSpecific(A,MAT_COOKIE);
139   PetscValidHeaderSpecific(B,MAT_COOKIE);
140   PetscCheckSameComm(ts,A);
141   PetscCheckSameComm(ts,B);
142   if (ts->problem_type == TS_NONLINEAR) {
143     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for nonlinear problems; use TSSetRHSJacobian()");
144   }
145 
146   ts->rhsmatrix = f;
147   ts->jacP      = ctx;
148   ts->A         = A;
149   ts->B         = B;
150 
151   PetscFunctionReturn(0);
152 }
153 
154 #undef __FUNC__
155 #define __FUNC__ "TSSetRHSJacobian"
156 /*@C
157    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
158    where U_t = F(U,t), as well as the location to store the matrix.
159 
160    Collective on TS
161 
162    Input Parameters:
163 +  ts  - the TS context obtained from TSCreate()
164 .  A   - Jacobian matrix
165 .  B   - preconditioner matrix (usually same as A)
166 .  f   - the Jacobian evaluation routine
167 -  ctx - [optional] user-defined context for private data for the
168          Jacobian evaluation routine (may be PETSC_NULL)
169 
170    Calling sequence of func:
171 $     func (TS ts,double t,Vec u,Mat *A,Mat *B,int *flag,void *ctx);
172 
173 +  t - current timestep
174 .  u - input vector
175 .  A - matrix A, where U_t = A(t)u
176 .  B - preconditioner matrix, usually the same as A
177 .  flag - flag indicating information about the preconditioner matrix
178           structure (same as flag in SLESSetOperators())
179 -  ctx - [optional] user-defined context for matrix evaluation routine
180 
181    Notes:
182    See SLESSetOperators() for important information about setting the flag
183    output parameter in the routine func().  Be sure to read this information!
184 
185    The routine func() takes Mat * as the matrix arguments rather than Mat.
186    This allows the matrix evaluation routine to replace A and/or B with a
187    completely new new matrix structure (not just different matrix elements)
188    when appropriate, for instance, if the nonzero structure is changing
189    throughout the global iterations.
190 
191    Level: beginner
192 
193 .keywords: TS, timestep, set, right-hand-side, Jacobian
194 
195 .seealso: TSDefaultComputeJacobianColor(),
196           SNESDefaultComputeJacobianColor()
197 
198 @*/
199 int TSSetRHSJacobian(TS ts,Mat A, Mat B,int (*f)(TS,double,Vec,Mat*,Mat*,
200                      MatStructure*,void*),void *ctx)
201 {
202   PetscFunctionBegin;
203   PetscValidHeaderSpecific(ts,TS_COOKIE);
204   PetscValidHeaderSpecific(A,MAT_COOKIE);
205   PetscValidHeaderSpecific(B,MAT_COOKIE);
206   PetscCheckSameComm(ts,A);
207   PetscCheckSameComm(ts,B);
208   if (ts->problem_type != TS_NONLINEAR) {
209     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for linear problems; use TSSetRHSMatrix()");
210   }
211 
212   ts->rhsjacobian = f;
213   ts->jacP        = ctx;
214   ts->A           = A;
215   ts->B           = B;
216   PetscFunctionReturn(0);
217 }
218 
219 #undef __FUNC__
220 #define __FUNC__ "TSComputeRHSBoundaryConditions"
221 /*
222    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
223 
224    Note: If the user did not provide a function but merely a matrix,
225    this routine applies the matrix.
226 */
227 int TSComputeRHSBoundaryConditions(TS ts,double t,Vec x)
228 {
229   int ierr;
230 
231   PetscFunctionBegin;
232   PetscValidHeaderSpecific(ts,TS_COOKIE);
233   PetscValidHeader(x);
234   PetscCheckSameComm(ts,x);
235 
236   if (ts->rhsbc) {
237     PetscStackPush("TS user boundary condition function");
238     ierr = (*ts->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
239     PetscStackPop;
240     PetscFunctionReturn(0);
241   }
242 
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNC__
247 #define __FUNC__ "TSSetRHSBoundaryConditions"
248 /*@C
249     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
250     boundary conditions for the function F.
251 
252     Collective on TS
253 
254     Input Parameters:
255 +   ts - the TS context obtained from TSCreate()
256 .   f - routine for evaluating the boundary condition function
257 -   ctx - [optional] user-defined context for private data for the
258           function evaluation routine (may be PETSC_NULL)
259 
260     Calling sequence of func:
261 $     func (TS ts,double t,Vec F,void *ctx);
262 
263 +   t - current timestep
264 .   F - function vector
265 -   ctx - [optional] user-defined function context
266 
267     Level: intermediate
268 
269 .keywords: TS, timestep, set, boundary conditions, function
270 @*/
271 int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,double,Vec,void*),void *ctx)
272 {
273   PetscFunctionBegin;
274 
275   PetscValidHeaderSpecific(ts,TS_COOKIE);
276   if (ts->problem_type != TS_LINEAR) {
277     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For linear problems only");
278   }
279   ts->rhsbc = f;
280   ts->bcP   = ctx;
281   PetscFunctionReturn(0);
282 }
283 
284 #undef __FUNC__
285 #define __FUNC__ "TSView"
286 /*@
287     TSView - Prints the TS data structure.
288 
289     Collective on TS, unless Viewer is VIEWER_STDOUT_SELF
290 
291     Input Parameters:
292 +   ts - the TS context obtained from TSCreate()
293 -   viewer - visualization context
294 
295     Options Database Key:
296 .   -ts_view - calls TSView() at end of TSStep()
297 
298     Notes:
299     The available visualization contexts include
300 +     VIEWER_STDOUT_SELF - standard output (default)
301 -     VIEWER_STDOUT_WORLD - synchronized standard
302          output where only the first processor opens
303          the file.  All other processors send their
304          data to the first processor to print.
305 
306     The user can open an alternative visualization context with
307     ViewerASCIIOpen() - output to a specified file.
308 
309     Level: beginner
310 
311 .keywords: TS, timestep, view
312 
313 .seealso: ViewerASCIIOpen()
314 @*/
315 int TSView(TS ts,Viewer viewer)
316 {
317   int        ierr;
318   char       *type;
319   PetscTruth isascii,isstring;
320 
321   PetscFunctionBegin;
322   PetscValidHeaderSpecific(ts,TS_COOKIE);
323   if (!viewer) viewer = VIEWER_STDOUT_SELF;
324   PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
325   PetscCheckSameComm(ts,viewer);
326 
327   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
328   ierr = PetscTypeCompare((PetscObject)viewer,STRING_VIEWER,&isstring);CHKERRQ(ierr);
329   if (isascii) {
330     ierr = ViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
331     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
332     if (type) {
333       ierr = ViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
334     } else {
335       ierr = ViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
336     }
337     if (ts->view) {
338       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
339       ierr = (*ts->view)(ts,viewer);CHKERRQ(ierr);
340       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
341     }
342     ierr = ViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
343     ierr = ViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
344     if (ts->problem_type == TS_NONLINEAR) {
345       ierr = ViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
346     }
347     ierr = ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
348   } else if (isstring) {
349     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
350     ierr = ViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
351   }
352   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
353   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
354   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
355   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
356   PetscFunctionReturn(0);
357 }
358 
359 
360 #undef __FUNC__
361 #define __FUNC__ "TSSetApplicationContext"
362 /*@C
363    TSSetApplicationContext - Sets an optional user-defined context for
364    the timesteppers.
365 
366    Collective on TS
367 
368    Input Parameters:
369 +  ts - the TS context obtained from TSCreate()
370 -  usrP - optional user context
371 
372    Level: intermediate
373 
374 .keywords: TS, timestep, set, application, context
375 
376 .seealso: TSGetApplicationContext()
377 @*/
378 int TSSetApplicationContext(TS ts,void *usrP)
379 {
380   PetscFunctionBegin;
381   PetscValidHeaderSpecific(ts,TS_COOKIE);
382   ts->user = usrP;
383   PetscFunctionReturn(0);
384 }
385 
386 #undef __FUNC__
387 #define __FUNC__ "TSGetApplicationContext"
388 /*@C
389     TSGetApplicationContext - Gets the user-defined context for the
390     timestepper.
391 
392     Not Collective
393 
394     Input Parameter:
395 .   ts - the TS context obtained from TSCreate()
396 
397     Output Parameter:
398 .   usrP - user context
399 
400     Level: intermediate
401 
402 .keywords: TS, timestep, get, application, context
403 
404 .seealso: TSSetApplicationContext()
405 @*/
406 int TSGetApplicationContext( TS ts,  void **usrP )
407 {
408   PetscFunctionBegin;
409   PetscValidHeaderSpecific(ts,TS_COOKIE);
410   *usrP = ts->user;
411   PetscFunctionReturn(0);
412 }
413 
414 #undef __FUNC__
415 #define __FUNC__ "TSGetTimeStepNumber"
416 /*@
417    TSGetTimeStepNumber - Gets the current number of timesteps.
418 
419    Not Collective
420 
421    Input Parameter:
422 .  ts - the TS context obtained from TSCreate()
423 
424    Output Parameter:
425 .  iter - number steps so far
426 
427    Level: intermediate
428 
429 .keywords: TS, timestep, get, iteration, number
430 @*/
431 int TSGetTimeStepNumber(TS ts,int* iter)
432 {
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(ts,TS_COOKIE);
435   *iter = ts->steps;
436   PetscFunctionReturn(0);
437 }
438 
439 #undef __FUNC__
440 #define __FUNC__ "TSSetInitialTimeStep"
441 /*@
442    TSSetInitialTimeStep - Sets the initial timestep to be used,
443    as well as the initial time.
444 
445    Collective on TS
446 
447    Input Parameters:
448 +  ts - the TS context obtained from TSCreate()
449 .  initial_time - the initial time
450 -  time_step - the size of the timestep
451 
452    Level: intermediate
453 
454 .seealso: TSSetTimeStep(), TSGetTimeStep()
455 
456 .keywords: TS, set, initial, timestep
457 @*/
458 int TSSetInitialTimeStep(TS ts,double initial_time,double time_step)
459 {
460   PetscFunctionBegin;
461   PetscValidHeaderSpecific(ts,TS_COOKIE);
462   ts->time_step         = time_step;
463   ts->initial_time_step = time_step;
464   ts->ptime             = initial_time;
465   PetscFunctionReturn(0);
466 }
467 
468 #undef __FUNC__
469 #define __FUNC__ "TSSetTimeStep"
470 /*@
471    TSSetTimeStep - Allows one to reset the timestep at any time,
472    useful for simple pseudo-timestepping codes.
473 
474    Collective on TS
475 
476    Input Parameters:
477 +  ts - the TS context obtained from TSCreate()
478 -  time_step - the size of the timestep
479 
480    Level: intermediate
481 
482 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
483 
484 .keywords: TS, set, timestep
485 @*/
486 int TSSetTimeStep(TS ts,double time_step)
487 {
488   PetscFunctionBegin;
489   PetscValidHeaderSpecific(ts,TS_COOKIE);
490   ts->time_step = time_step;
491   PetscFunctionReturn(0);
492 }
493 
494 #undef __FUNC__
495 #define __FUNC__ "TSGetTimeStep"
496 /*@
497    TSGetTimeStep - Gets the current timestep size.
498 
499    Not Collective
500 
501    Input Parameter:
502 .  ts - the TS context obtained from TSCreate()
503 
504    Output Parameter:
505 .  dt - the current timestep size
506 
507    Level: intermediate
508 
509 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
510 
511 .keywords: TS, get, timestep
512 @*/
513 int TSGetTimeStep(TS ts,double* dt)
514 {
515   PetscFunctionBegin;
516   PetscValidHeaderSpecific(ts,TS_COOKIE);
517   *dt = ts->time_step;
518   PetscFunctionReturn(0);
519 }
520 
521 #undef __FUNC__
522 #define __FUNC__ "TSGetSolution"
523 /*@C
524    TSGetSolution - Returns the solution at the present timestep. It
525    is valid to call this routine inside the function that you are evaluating
526    in order to move to the new timestep. This vector not changed until
527    the solution at the next timestep has been calculated.
528 
529    Not Collective, but Vec returned is parallel if TS is parallel
530 
531    Input Parameter:
532 .  ts - the TS context obtained from TSCreate()
533 
534    Output Parameter:
535 .  v - the vector containing the solution
536 
537    Level: intermediate
538 
539 .seealso: TSGetTimeStep()
540 
541 .keywords: TS, timestep, get, solution
542 @*/
543 int TSGetSolution(TS ts,Vec *v)
544 {
545   PetscFunctionBegin;
546   PetscValidHeaderSpecific(ts,TS_COOKIE);
547   *v = ts->vec_sol_always;
548   PetscFunctionReturn(0);
549 }
550 
551 #undef __FUNC__
552 #define __FUNC__ "TSPublish_Petsc"
553 static int TSPublish_Petsc(PetscObject obj)
554 {
555 #if defined(PETSC_HAVE_AMS)
556   TS   v = (TS) obj;
557   int  ierr;
558 #endif
559 
560   PetscFunctionBegin;
561 
562 #if defined(PETSC_HAVE_AMS)
563   /* if it is already published then return */
564   if (v->amem >=0 ) PetscFunctionReturn(0);
565 
566   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
567   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ,
568                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
569   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ,
570                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
571   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1,
572                                AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
573   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
574 #endif
575   PetscFunctionReturn(0);
576 }
577 
578 /* -----------------------------------------------------------*/
579 
580 #undef __FUNC__
581 #define __FUNC__ "TSCreate"
582 /*@C
583    TSCreate - Creates a timestepper context.
584 
585    Collective on MPI_Comm
586 
587    Input Parameter:
588 +  comm - MPI communicator
589 -  type - One of  TS_LINEAR,TS_NONLINEAR
590    where these types refer to problems of the forms
591 .vb
592          U_t = A U
593          U_t = A(t) U
594          U_t = F(t,U)
595 .ve
596 
597    Output Parameter:
598 .  outts - the new TS context
599 
600    Level: beginner
601 
602 .keywords: TS, timestep, create, context
603 
604 .seealso: TSSetUp(), TSStep(), TSDestroy()
605 @*/
606 int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts)
607 {
608   TS   ts;
609 
610   PetscFunctionBegin;
611   *outts = 0;
612   PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView);
613   PLogObjectCreate(ts);
614   ts->bops->publish     = TSPublish_Petsc;
615   ts->max_steps         = 5000;
616   ts->max_time          = 5.0;
617   ts->time_step         = .1;
618   ts->initial_time_step = ts->time_step;
619   ts->steps             = 0;
620   ts->ptime             = 0.0;
621   ts->data              = 0;
622   ts->view              = 0;
623   ts->setupcalled       = 0;
624   ts->problem_type      = problemtype;
625   ts->numbermonitors    = 0;
626   ts->linear_its        = 0;
627   ts->nonlinear_its     = 0;
628 
629   *outts = ts;
630   PetscFunctionReturn(0);
631 }
632 
633 /* ----- Routines to initialize and destroy a timestepper ---- */
634 
635 #undef __FUNC__
636 #define __FUNC__ "TSSetUp"
637 /*@
638    TSSetUp - Sets up the internal data structures for the later use
639    of a timestepper.
640 
641    Collective on TS
642 
643    Input Parameter:
644 .  ts - the TS context obtained from TSCreate()
645 
646    Notes:
647    For basic use of the TS solvers the user need not explicitly call
648    TSSetUp(), since these actions will automatically occur during
649    the call to TSStep().  However, if one wishes to control this
650    phase separately, TSSetUp() should be called after TSCreate()
651    and optional routines of the form TSSetXXX(), but before TSStep().
652 
653    Level: advanced
654 
655 .keywords: TS, timestep, setup
656 
657 .seealso: TSCreate(), TSStep(), TSDestroy()
658 @*/
659 int TSSetUp(TS ts)
660 {
661   int ierr;
662 
663   PetscFunctionBegin;
664   PetscValidHeaderSpecific(ts,TS_COOKIE);
665   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call TSSetSolution() first");
666   if (!ts->type_name) {
667     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
668   }
669   ierr = (*ts->setup)(ts);CHKERRQ(ierr);
670   ts->setupcalled = 1;
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNC__
675 #define __FUNC__ "TSDestroy"
676 /*@C
677    TSDestroy - Destroys the timestepper context that was created
678    with TSCreate().
679 
680    Collective on TS
681 
682    Input Parameter:
683 .  ts - the TS context obtained from TSCreate()
684 
685    Level: beginner
686 
687 .keywords: TS, timestepper, destroy
688 
689 .seealso: TSCreate(), TSSetUp(), TSSolve()
690 @*/
691 int TSDestroy(TS ts)
692 {
693   int ierr;
694 
695   PetscFunctionBegin;
696   PetscValidHeaderSpecific(ts,TS_COOKIE);
697   if (--ts->refct > 0) PetscFunctionReturn(0);
698 
699   /* if memory was published with AMS then destroy it */
700   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
701 
702   if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);}
703   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
704   ierr = (*(ts)->destroy)(ts);CHKERRQ(ierr);
705   PLogObjectDestroy((PetscObject)ts);
706   PetscHeaderDestroy((PetscObject)ts);
707   PetscFunctionReturn(0);
708 }
709 
710 #undef __FUNC__
711 #define __FUNC__ "TSGetSNES"
712 /*@C
713    TSGetSNES - Returns the SNES (nonlinear solver) associated with
714    a TS (timestepper) context. Valid only for nonlinear problems.
715 
716    Not Collective, but SNES is parallel if TS is parallel
717 
718    Input Parameter:
719 .  ts - the TS context obtained from TSCreate()
720 
721    Output Parameter:
722 .  snes - the nonlinear solver context
723 
724    Notes:
725    The user can then directly manipulate the SNES context to set various
726    options, etc.  Likewise, the user can then extract and manipulate the
727    SLES, KSP, and PC contexts as well.
728 
729    TSGetSNES() does not work for integrators that do not use SNES; in
730    this case TSGetSNES() returns PETSC_NULL in snes.
731 
732    Level: beginner
733 
734 .keywords: timestep, get, SNES
735 @*/
736 int TSGetSNES(TS ts,SNES *snes)
737 {
738   PetscFunctionBegin;
739   PetscValidHeaderSpecific(ts,TS_COOKIE);
740   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Nonlinear only; use TSGetSLES()");
741   *snes = ts->snes;
742   PetscFunctionReturn(0);
743 }
744 
745 #undef __FUNC__
746 #define __FUNC__ "TSGetSLES"
747 /*@C
748    TSGetSLES - Returns the SLES (linear solver) associated with
749    a TS (timestepper) context.
750 
751    Not Collective, but SLES is parallel if TS is parallel
752 
753    Input Parameter:
754 .  ts - the TS context obtained from TSCreate()
755 
756    Output Parameter:
757 .  sles - the nonlinear solver context
758 
759    Notes:
760    The user can then directly manipulate the SLES context to set various
761    options, etc.  Likewise, the user can then extract and manipulate the
762    KSP and PC contexts as well.
763 
764    TSGetSLES() does not work for integrators that do not use SLES;
765    in this case TSGetSLES() returns PETSC_NULL in sles.
766 
767    Level: beginner
768 
769 .keywords: timestep, get, SLES
770 @*/
771 int TSGetSLES(TS ts,SLES *sles)
772 {
773   PetscFunctionBegin;
774   PetscValidHeaderSpecific(ts,TS_COOKIE);
775   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Linear only; use TSGetSNES()");
776   *sles = ts->sles;
777   PetscFunctionReturn(0);
778 }
779 
780 /* ----------- Routines to set solver parameters ---------- */
781 
782 #undef __FUNC__
783 #define __FUNC__ "TSSetDuration"
784 /*@
785    TSSetDuration - Sets the maximum number of timesteps to use and
786    maximum time for iteration.
787 
788    Collective on TS
789 
790    Input Parameters:
791 +  ts - the TS context obtained from TSCreate()
792 .  maxsteps - maximum number of iterations to use
793 -  maxtime - final time to iterate to
794 
795    Options Database Keys:
796 .  -ts_max_steps <maxsteps> - Sets maxsteps
797 .  -ts_max_time <maxtime> - Sets maxtime
798 
799    Notes:
800    The default maximum number of iterations is 5000. Default time is 5.0
801 
802    Level: intermediate
803 
804 .keywords: TS, timestep, set, maximum, iterations
805 @*/
806 int TSSetDuration(TS ts,int maxsteps,double maxtime)
807 {
808   PetscFunctionBegin;
809   PetscValidHeaderSpecific(ts,TS_COOKIE);
810   ts->max_steps = maxsteps;
811   ts->max_time  = maxtime;
812   PetscFunctionReturn(0);
813 }
814 
815 #undef __FUNC__
816 #define __FUNC__ "TSSetSolution"
817 /*@
818    TSSetSolution - Sets the initial solution vector
819    for use by the TS routines.
820 
821    Collective on TS and Vec
822 
823    Input Parameters:
824 +  ts - the TS context obtained from TSCreate()
825 -  x - the solution vector
826 
827    Level: beginner
828 
829 .keywords: TS, timestep, set, solution, initial conditions
830 @*/
831 int TSSetSolution(TS ts,Vec x)
832 {
833   PetscFunctionBegin;
834   PetscValidHeaderSpecific(ts,TS_COOKIE);
835   ts->vec_sol        = ts->vec_sol_always = x;
836   PetscFunctionReturn(0);
837 }
838 
839 /* ------------ Routines to set performance monitoring options ----------- */
840 
841 #undef __FUNC__
842 #define __FUNC__ "TSSetMonitor"
843 /*@C
844    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
845    timestep to display the iteration's  progress.
846 
847    Collective on TS
848 
849    Input Parameters:
850 +  ts - the TS context obtained from TSCreate()
851 .  func - monitoring routine
852 -  mctx - [optional] user-defined context for private data for the
853           monitor routine (may be PETSC_NULL)
854 
855    Calling sequence of func:
856 $    int func(TS ts,int steps,double time,Vec x,void *mctx)
857 
858 +    ts - the TS context
859 .    steps - iteration number
860 .    time - current timestep
861 .    x - current iterate
862 -    mctx - [optional] monitoring context
863 
864    Notes:
865    This routine adds an additional monitor to the list of monitors that
866    already has been loaded.
867 
868    Level: intermediate
869 
870 .keywords: TS, timestep, set, monitor
871 
872 .seealso: TSDefaultMonitor(), TSClearMonitor()
873 @*/
874 int TSSetMonitor(TS ts, int (*monitor)(TS,int,double,Vec,void*), void *mctx )
875 {
876   PetscFunctionBegin;
877   PetscValidHeaderSpecific(ts,TS_COOKIE);
878   if (ts->numbermonitors >= MAXTSMONITORS) {
879     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
880   }
881   ts->monitor[ts->numbermonitors]           = monitor;
882   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
883   PetscFunctionReturn(0);
884 }
885 
886 #undef __FUNC__
887 #define __FUNC__ "TSClearMonitor"
888 /*@C
889    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
890 
891    Collective on TS
892 
893    Input Parameters:
894 .  ts - the TS context obtained from TSCreate()
895 
896    Notes:
897    There is no way to remove a single, specific monitor.
898 
899    Level: intermediate
900 
901 .keywords: TS, timestep, set, monitor
902 
903 .seealso: TSDefaultMonitor(), TSSetMonitor()
904 @*/
905 int TSClearMonitor(TS ts)
906 {
907   PetscFunctionBegin;
908   PetscValidHeaderSpecific(ts,TS_COOKIE);
909   ts->numbermonitors = 0;
910   PetscFunctionReturn(0);
911 }
912 
913 #undef __FUNC__
914 #define __FUNC__ "TSDefaultMonitor"
915 int TSDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
916 {
917   int ierr;
918 
919   PetscFunctionBegin;
920   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,time);CHKERRQ(ierr);
921   PetscFunctionReturn(0);
922 }
923 
924 #undef __FUNC__
925 #define __FUNC__ "TSStep"
926 /*@
927    TSStep - Steps the requested number of timesteps.
928 
929    Collective on TS
930 
931    Input Parameter:
932 .  ts - the TS context obtained from TSCreate()
933 
934    Output Parameters:
935 +  steps - number of iterations until termination
936 -  time - time until termination
937 
938    Level: beginner
939 
940 .keywords: TS, timestep, solve
941 
942 .seealso: TSCreate(), TSSetUp(), TSDestroy()
943 @*/
944 int TSStep(TS ts,int *steps,double *time)
945 {
946   int ierr,flg;
947 
948   PetscFunctionBegin;
949   PetscValidHeaderSpecific(ts,TS_COOKIE);
950   if (!ts->setupcalled) {ierr = TSSetUp(ts);CHKERRQ(ierr);}
951   PLogEventBegin(TS_Step,ts,0,0,0);
952   ierr = (*(ts)->step)(ts,steps,time);CHKERRQ(ierr);
953   PLogEventEnd(TS_Step,ts,0,0,0);
954   ierr = OptionsHasName(ts->prefix,"-ts_view",&flg);CHKERRQ(ierr);
955   if (flg) {
956     ierr = TSView(ts,VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
957   }
958   PetscFunctionReturn(0);
959 }
960 
961 #undef __FUNC__
962 #define __FUNC__ "TSMonitor"
963 /*
964      Runs the user provided monitor routines, if they exists.
965 */
966 int TSMonitor(TS ts,int step,double time,Vec x)
967 {
968   int i,ierr,n = ts->numbermonitors;
969 
970   PetscFunctionBegin;
971   for ( i=0; i<n; i++ ) {
972     ierr = (*ts->monitor[i])(ts,step,time,x,ts->monitorcontext[i]);CHKERRQ(ierr);
973   }
974   PetscFunctionReturn(0);
975 }
976 
977 /* ------------------------------------------------------------------------*/
978 
979 #undef __FUNC__
980 #define __FUNC__ "TSLGMonitorCreate"
981 /*@C
982    TSLGMonitorCreate - Creates a line graph context for use with
983    TS to monitor convergence of preconditioned residual norms.
984 
985    Collective on TS
986 
987    Input Parameters:
988 +  host - the X display to open, or null for the local machine
989 .  label - the title to put in the title bar
990 .  x, y - the screen coordinates of the upper left coordinate of
991           the window
992 -  m, n - the screen width and height in pixels
993 
994    Output Parameter:
995 .  draw - the drawing context
996 
997    Options Database Key:
998 .  -ts_xmonitor - automatically sets line graph monitor
999 
1000    Notes:
1001    Use TSLGMonitorDestroy() to destroy this line graph, not DrawLGDestroy().
1002 
1003    Level: intermediate
1004 
1005 .keywords: TS, monitor, line graph, residual, create
1006 
1007 .seealso: TSLGMonitorDestroy(), TSSetMonitor()
1008 @*/
1009 int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,
1010                        int n, DrawLG *draw)
1011 {
1012   Draw win;
1013   int  ierr;
1014 
1015   PetscFunctionBegin;
1016   ierr = DrawOpenX(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1017   ierr = DrawLGCreate(win,1,draw);CHKERRQ(ierr);
1018   ierr = DrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1019 
1020   PLogObjectParent(*draw,win);
1021   PetscFunctionReturn(0);
1022 }
1023 
1024 #undef __FUNC__
1025 #define __FUNC__ "TSLGMonitor"
1026 int TSLGMonitor(TS ts,int n,double time,Vec v,void *monctx)
1027 {
1028   DrawLG lg = (DrawLG) monctx;
1029   double x,y = time;
1030   int    ierr;
1031 
1032   PetscFunctionBegin;
1033   if (!n) {ierr = DrawLGReset(lg);CHKERRQ(ierr);}
1034   x = (double) n;
1035   ierr = DrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1036   if (n < 20 || (n % 5)) {
1037     ierr = DrawLGDraw(lg);CHKERRQ(ierr);
1038   }
1039   PetscFunctionReturn(0);
1040 }
1041 
1042 #undef __FUNC__
1043 #define __FUNC__ "TSLGMonitorDestroy"
1044 /*@C
1045    TSLGMonitorDestroy - Destroys a line graph context that was created
1046    with TSLGMonitorCreate().
1047 
1048    Collective on DrawLG
1049 
1050    Input Parameter:
1051 .  draw - the drawing context
1052 
1053    Level: intermediate
1054 
1055 .keywords: TS, monitor, line graph, destroy
1056 
1057 .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1058 @*/
1059 int TSLGMonitorDestroy(DrawLG drawlg)
1060 {
1061   Draw draw;
1062   int  ierr;
1063 
1064   PetscFunctionBegin;
1065   ierr = DrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1066   ierr = DrawDestroy(draw);CHKERRQ(ierr);
1067   ierr = DrawLGDestroy(drawlg);CHKERRQ(ierr);
1068   PetscFunctionReturn(0);
1069 }
1070 
1071 #undef __FUNC__
1072 #define __FUNC__ "TSGetTime"
1073 /*@
1074    TSGetTime - Gets the current time.
1075 
1076    Not Collective
1077 
1078    Input Parameter:
1079 .  ts - the TS context obtained from TSCreate()
1080 
1081    Output Parameter:
1082 .  t  - the current time
1083 
1084    Contributed by: Matthew Knepley
1085 
1086    Level: beginner
1087 
1088 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1089 
1090 .keywords: TS, get, time
1091 @*/
1092 int TSGetTime(TS ts, double* t)
1093 {
1094   PetscFunctionBegin;
1095   PetscValidHeaderSpecific(ts, TS_COOKIE);
1096   *t = ts->ptime;
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 #undef __FUNC__
1101 #define __FUNC__ "TSGetProblemType"
1102 /*@C
1103    TSGetProblemType - Returns the problem type of a TS (timestepper) context.
1104 
1105    Not Collective
1106 
1107    Input Parameter:
1108 .  ts   - The TS context obtained from TSCreate()
1109 
1110    Output Parameter:
1111 .  type - The problem type, TS_LINEAR or TS_NONLINEAR
1112 
1113    Level: intermediate
1114 
1115    Contributed by: Matthew Knepley
1116 
1117 .keywords: ts, get, type
1118 
1119 @*/
1120 int TSGetProblemType(TS ts, TSProblemType *type)
1121 {
1122   PetscFunctionBegin;
1123   PetscValidHeaderSpecific(ts, TS_COOKIE);
1124   *type = ts->problem_type;
1125   PetscFunctionReturn(0);
1126 }
1127 
1128 #undef __FUNC__
1129 #define __FUNC__ "TSSetOptionsPrefix"
1130 /*@C
1131    TSSetOptionsPrefix - Sets the prefix used for searching for all
1132    TS options in the database.
1133 
1134    Collective on TS
1135 
1136    Input Parameter:
1137 +  ts     - The TS context
1138 -  prefix - The prefix to prepend to all option names
1139 
1140    Notes:
1141    A hyphen (-) must NOT be given at the beginning of the prefix name.
1142    The first character of all runtime options is AUTOMATICALLY the
1143    hyphen.
1144 
1145    Contributed by: Matthew Knepley
1146 
1147    Level: advanced
1148 
1149 .keywords: TS, set, options, prefix, database
1150 
1151 .seealso: TSSetFromOptions()
1152 
1153 @*/
1154 int TSSetOptionsPrefix(TS ts, char *prefix)
1155 {
1156   int ierr;
1157 
1158   PetscFunctionBegin;
1159   PetscValidHeaderSpecific(ts, TS_COOKIE);
1160   ierr = PetscObjectSetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1161   switch(ts->problem_type) {
1162     case TS_NONLINEAR:
1163       ierr = SNESSetOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr);
1164       break;
1165     case TS_LINEAR:
1166       ierr = SLESSetOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr);
1167       break;
1168   }
1169   PetscFunctionReturn(0);
1170 }
1171 
1172 
1173 #undef __FUNC__
1174 #define __FUNC__ "TSAppendOptionsPrefix"
1175 /*@C
1176    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1177    TS options in the database.
1178 
1179    Collective on TS
1180 
1181    Input Parameter:
1182 +  ts     - The TS context
1183 -  prefix - The prefix to prepend to all option names
1184 
1185    Notes:
1186    A hyphen (-) must NOT be given at the beginning of the prefix name.
1187    The first character of all runtime options is AUTOMATICALLY the
1188    hyphen.
1189 
1190    Contributed by: Matthew Knepley
1191 
1192    Level: advanced
1193 
1194 .keywords: TS, append, options, prefix, database
1195 
1196 .seealso: TSGetOptionsPrefix()
1197 
1198 @*/
1199 int TSAppendOptionsPrefix(TS ts, char *prefix)
1200 {
1201   int ierr;
1202 
1203   PetscFunctionBegin;
1204   PetscValidHeaderSpecific(ts, TS_COOKIE);
1205   ierr = PetscObjectAppendOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1206   switch(ts->problem_type) {
1207     case TS_NONLINEAR:
1208       ierr = SNESAppendOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr);
1209       break;
1210     case TS_LINEAR:
1211       ierr = SLESAppendOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr);
1212       break;
1213   }
1214   PetscFunctionReturn(0);
1215 }
1216 
1217 #undef __FUNC__
1218 #define __FUNC__ "TSGetOptionsPrefix"
1219 /*@C
1220    TSGetOptionsPrefix - Sets the prefix used for searching for all
1221    TS options in the database.
1222 
1223    Not Collective
1224 
1225    Input Parameter:
1226 .  ts - The TS context
1227 
1228    Output Parameter:
1229 .  prefix - A pointer to the prefix string used
1230 
1231    Contributed by: Matthew Knepley
1232 
1233    Notes: On the fortran side, the user should pass in a string 'prifix' of
1234    sufficient length to hold the prefix.
1235 
1236    Level: intermediate
1237 
1238 .keywords: TS, get, options, prefix, database
1239 
1240 .seealso: TSAppendOptionsPrefix()
1241 @*/
1242 int TSGetOptionsPrefix(TS ts, char **prefix)
1243 {
1244   int ierr;
1245 
1246   PetscFunctionBegin;
1247   PetscValidHeaderSpecific(ts, TS_COOKIE);
1248   ierr = PetscObjectGetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1249   PetscFunctionReturn(0);
1250 }
1251 
1252 #undef __FUNC__
1253 #define __FUNC__ "TSGetRHSMatrix"
1254 /*@C
1255    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1256 
1257    Not Collective, but parallel objects are returned if TS is parallel
1258 
1259    Input Parameter:
1260 .  ts  - The TS context obtained from TSCreate()
1261 
1262    Output Parameters:
1263 +  A   - The matrix A, where U_t = A(t) U
1264 .  M   - The preconditioner matrix, usually the same as A
1265 -  ctx - User-defined context for matrix evaluation routine
1266 
1267    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1268 
1269    Contributed by: Matthew Knepley
1270 
1271    Level: intermediate
1272 
1273 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1274 
1275 .keywords: TS, timestep, get, matrix
1276 
1277 @*/
1278 int TSGetRHSMatrix(TS ts, Mat *A, Mat *M, void **ctx)
1279 {
1280   PetscFunctionBegin;
1281   PetscValidHeaderSpecific(ts, TS_COOKIE);
1282   if (A)   *A = ts->A;
1283   if (M)   *M = ts->B;
1284   if (ctx) *ctx = ts->jacP;
1285   PetscFunctionReturn(0);
1286 }
1287 
1288 #undef __FUNC__
1289 #define __FUNC__ "TSGetRHSJacobian"
1290 /*@C
1291    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1292 
1293    Not Collective, but parallel objects are returned if TS is parallel
1294 
1295    Input Parameter:
1296 .  ts  - The TS context obtained from TSCreate()
1297 
1298    Output Parameters:
1299 +  J   - The Jacobian J of F, where U_t = F(U,t)
1300 .  M   - The preconditioner matrix, usually the same as J
1301 - ctx - User-defined context for Jacobian evaluation routine
1302 
1303    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1304 
1305    Contributed by: Matthew Knepley
1306 
1307    Level: intermediate
1308 
1309 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1310 
1311 .keywords: TS, timestep, get, matrix, Jacobian
1312 @*/
1313 int TSGetRHSJacobian(TS ts, Mat *J, Mat *M, void **ctx)
1314 {
1315   int ierr;
1316 
1317   PetscFunctionBegin;
1318   ierr = TSGetRHSMatrix(ts, J, M, ctx);CHKERRQ(ierr);
1319   PetscFunctionReturn(0);
1320 }
1321 
1322 /*MC
1323    TSRegister - Adds a method to the timestepping solver package.
1324 
1325    Synopsis:
1326 
1327    TSRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(TS))
1328 
1329    Not collective
1330 
1331    Input Parameters:
1332 +  name_solver - name of a new user-defined solver
1333 .  path - path (either absolute or relative) the library containing this solver
1334 .  name_create - name of routine to create method context
1335 -  routine_create - routine to create method context
1336 
1337    Notes:
1338    TSRegister() may be called multiple times to add several user-defined solvers.
1339 
1340    If dynamic libraries are used, then the fourth input argument (routine_create)
1341    is ignored.
1342 
1343    Sample usage:
1344 .vb
1345    TSRegister("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
1346               "MySolverCreate",MySolverCreate);
1347 .ve
1348 
1349    Then, your solver can be chosen with the procedural interface via
1350 $     TSSetType(ts,"my_solver")
1351    or at runtime via the option
1352 $     -ts_type my_solver
1353 
1354    Level: advanced
1355 
1356    $PETSC_ARCH, $PETSC_DIR, $PETSC_LDIR, and $BOPT occuring in pathname will be replaced with appropriate values.
1357 
1358 .keywords: TS, register
1359 
1360 .seealso: TSRegisterAll(), TSRegisterDestroy()
1361 M*/
1362 
1363 #undef __FUNC__
1364 #define __FUNC__ "TSRegister_Private"
1365 int TSRegister_Private(char *sname,char *path,char *name,int (*function)(TS))
1366 {
1367   char fullname[256];
1368   int  ierr;
1369 
1370   PetscFunctionBegin;
1371   ierr = FListConcat_Private(path,name,fullname); CHKERRQ(ierr);
1372   FListAdd_Private(&TSList,sname,fullname,        (int (*)(void*))function);
1373   PetscFunctionReturn(0);
1374 }
1375