xref: /petsc/src/ts/interface/ts.c (revision 0f5bd95cca42d693ada6d048329ab39533680180)
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   int      isascii,isstring;
320 
321   PetscFunctionBegin;
322   PetscValidHeaderSpecific(ts,TS_COOKIE);
323   if (!viewer) viewer = VIEWER_STDOUT_SELF;
324   PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
325 
326   isascii = PetscTypeCompare(viewer,ASCII_VIEWER);
327   isstring = PetscTypeCompare(viewer,STRING_VIEWER);
328   if (isascii) {
329     ierr = ViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
330     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
331     if (type) {
332       ierr = ViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
333     } else {
334       ierr = ViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
335     }
336     if (ts->view) {
337       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
338       ierr = (*ts->view)(ts,viewer);CHKERRQ(ierr);
339       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
340     }
341     ierr = ViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
342     ierr = ViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
343     if (ts->problem_type == TS_NONLINEAR) {
344       ierr = ViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
345     }
346     ierr = ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
347   } else if (isstring) {
348     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
349     ierr = ViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
350   }
351   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
352   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
353   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
354   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
355   PetscFunctionReturn(0);
356 }
357 
358 
359 #undef __FUNC__
360 #define __FUNC__ "TSSetApplicationContext"
361 /*@C
362    TSSetApplicationContext - Sets an optional user-defined context for
363    the timesteppers.
364 
365    Collective on TS
366 
367    Input Parameters:
368 +  ts - the TS context obtained from TSCreate()
369 -  usrP - optional user context
370 
371    Level: intermediate
372 
373 .keywords: TS, timestep, set, application, context
374 
375 .seealso: TSGetApplicationContext()
376 @*/
377 int TSSetApplicationContext(TS ts,void *usrP)
378 {
379   PetscFunctionBegin;
380   PetscValidHeaderSpecific(ts,TS_COOKIE);
381   ts->user = usrP;
382   PetscFunctionReturn(0);
383 }
384 
385 #undef __FUNC__
386 #define __FUNC__ "TSGetApplicationContext"
387 /*@C
388     TSGetApplicationContext - Gets the user-defined context for the
389     timestepper.
390 
391     Not Collective
392 
393     Input Parameter:
394 .   ts - the TS context obtained from TSCreate()
395 
396     Output Parameter:
397 .   usrP - user context
398 
399     Level: intermediate
400 
401 .keywords: TS, timestep, get, application, context
402 
403 .seealso: TSSetApplicationContext()
404 @*/
405 int TSGetApplicationContext( TS ts,  void **usrP )
406 {
407   PetscFunctionBegin;
408   PetscValidHeaderSpecific(ts,TS_COOKIE);
409   *usrP = ts->user;
410   PetscFunctionReturn(0);
411 }
412 
413 #undef __FUNC__
414 #define __FUNC__ "TSGetTimeStepNumber"
415 /*@
416    TSGetTimeStepNumber - Gets the current number of timesteps.
417 
418    Not Collective
419 
420    Input Parameter:
421 .  ts - the TS context obtained from TSCreate()
422 
423    Output Parameter:
424 .  iter - number steps so far
425 
426    Level: intermediate
427 
428 .keywords: TS, timestep, get, iteration, number
429 @*/
430 int TSGetTimeStepNumber(TS ts,int* iter)
431 {
432   PetscFunctionBegin;
433   PetscValidHeaderSpecific(ts,TS_COOKIE);
434   *iter = ts->steps;
435   PetscFunctionReturn(0);
436 }
437 
438 #undef __FUNC__
439 #define __FUNC__ "TSSetInitialTimeStep"
440 /*@
441    TSSetInitialTimeStep - Sets the initial timestep to be used,
442    as well as the initial time.
443 
444    Collective on TS
445 
446    Input Parameters:
447 +  ts - the TS context obtained from TSCreate()
448 .  initial_time - the initial time
449 -  time_step - the size of the timestep
450 
451    Level: intermediate
452 
453 .seealso: TSSetTimeStep(), TSGetTimeStep()
454 
455 .keywords: TS, set, initial, timestep
456 @*/
457 int TSSetInitialTimeStep(TS ts,double initial_time,double time_step)
458 {
459   PetscFunctionBegin;
460   PetscValidHeaderSpecific(ts,TS_COOKIE);
461   ts->time_step         = time_step;
462   ts->initial_time_step = time_step;
463   ts->ptime             = initial_time;
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNC__
468 #define __FUNC__ "TSSetTimeStep"
469 /*@
470    TSSetTimeStep - Allows one to reset the timestep at any time,
471    useful for simple pseudo-timestepping codes.
472 
473    Collective on TS
474 
475    Input Parameters:
476 +  ts - the TS context obtained from TSCreate()
477 -  time_step - the size of the timestep
478 
479    Level: intermediate
480 
481 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
482 
483 .keywords: TS, set, timestep
484 @*/
485 int TSSetTimeStep(TS ts,double time_step)
486 {
487   PetscFunctionBegin;
488   PetscValidHeaderSpecific(ts,TS_COOKIE);
489   ts->time_step = time_step;
490   PetscFunctionReturn(0);
491 }
492 
493 #undef __FUNC__
494 #define __FUNC__ "TSGetTimeStep"
495 /*@
496    TSGetTimeStep - Gets the current timestep size.
497 
498    Not Collective
499 
500    Input Parameter:
501 .  ts - the TS context obtained from TSCreate()
502 
503    Output Parameter:
504 .  dt - the current timestep size
505 
506    Level: intermediate
507 
508 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
509 
510 .keywords: TS, get, timestep
511 @*/
512 int TSGetTimeStep(TS ts,double* dt)
513 {
514   PetscFunctionBegin;
515   PetscValidHeaderSpecific(ts,TS_COOKIE);
516   *dt = ts->time_step;
517   PetscFunctionReturn(0);
518 }
519 
520 #undef __FUNC__
521 #define __FUNC__ "TSGetSolution"
522 /*@C
523    TSGetSolution - Returns the solution at the present timestep. It
524    is valid to call this routine inside the function that you are evaluating
525    in order to move to the new timestep. This vector not changed until
526    the solution at the next timestep has been calculated.
527 
528    Not Collective, but Vec returned is parallel if TS is parallel
529 
530    Input Parameter:
531 .  ts - the TS context obtained from TSCreate()
532 
533    Output Parameter:
534 .  v - the vector containing the solution
535 
536    Level: intermediate
537 
538 .seealso: TSGetTimeStep()
539 
540 .keywords: TS, timestep, get, solution
541 @*/
542 int TSGetSolution(TS ts,Vec *v)
543 {
544   PetscFunctionBegin;
545   PetscValidHeaderSpecific(ts,TS_COOKIE);
546   *v = ts->vec_sol_always;
547   PetscFunctionReturn(0);
548 }
549 
550 #undef __FUNC__
551 #define __FUNC__ "TSPublish_Petsc"
552 static int TSPublish_Petsc(PetscObject obj)
553 {
554 #if defined(PETSC_HAVE_AMS)
555   TS   v = (TS) obj;
556   int  ierr;
557 #endif
558 
559   PetscFunctionBegin;
560 
561 #if defined(PETSC_HAVE_AMS)
562   /* if it is already published then return */
563   if (v->amem >=0 ) PetscFunctionReturn(0);
564 
565   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
566   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ,
567                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
568   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ,
569                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
570   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1,
571                                AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
572   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
573 #endif
574   PetscFunctionReturn(0);
575 }
576 
577 /* -----------------------------------------------------------*/
578 
579 #undef __FUNC__
580 #define __FUNC__ "TSCreate"
581 /*@C
582    TSCreate - Creates a timestepper context.
583 
584    Collective on MPI_Comm
585 
586    Input Parameter:
587 +  comm - MPI communicator
588 -  type - One of  TS_LINEAR,TS_NONLINEAR
589    where these types refer to problems of the forms
590 .vb
591          U_t = A U
592          U_t = A(t) U
593          U_t = F(t,U)
594 .ve
595 
596    Output Parameter:
597 .  outts - the new TS context
598 
599    Level: beginner
600 
601 .keywords: TS, timestep, create, context
602 
603 .seealso: TSSetUp(), TSStep(), TSDestroy()
604 @*/
605 int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts)
606 {
607   TS   ts;
608 
609   PetscFunctionBegin;
610   *outts = 0;
611   PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView);
612   PLogObjectCreate(ts);
613   ts->bops->publish     = TSPublish_Petsc;
614   ts->max_steps         = 5000;
615   ts->max_time          = 5.0;
616   ts->time_step         = .1;
617   ts->initial_time_step = ts->time_step;
618   ts->steps             = 0;
619   ts->ptime             = 0.0;
620   ts->data              = 0;
621   ts->view              = 0;
622   ts->setupcalled       = 0;
623   ts->problem_type      = problemtype;
624   ts->numbermonitors    = 0;
625   ts->linear_its        = 0;
626   ts->nonlinear_its     = 0;
627 
628   *outts = ts;
629   PetscFunctionReturn(0);
630 }
631 
632 /* ----- Routines to initialize and destroy a timestepper ---- */
633 
634 #undef __FUNC__
635 #define __FUNC__ "TSSetUp"
636 /*@
637    TSSetUp - Sets up the internal data structures for the later use
638    of a timestepper.
639 
640    Collective on TS
641 
642    Input Parameter:
643 .  ts - the TS context obtained from TSCreate()
644 
645    Notes:
646    For basic use of the TS solvers the user need not explicitly call
647    TSSetUp(), since these actions will automatically occur during
648    the call to TSStep().  However, if one wishes to control this
649    phase separately, TSSetUp() should be called after TSCreate()
650    and optional routines of the form TSSetXXX(), but before TSStep().
651 
652    Level: advanced
653 
654 .keywords: TS, timestep, setup
655 
656 .seealso: TSCreate(), TSStep(), TSDestroy()
657 @*/
658 int TSSetUp(TS ts)
659 {
660   int ierr;
661 
662   PetscFunctionBegin;
663   PetscValidHeaderSpecific(ts,TS_COOKIE);
664   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call TSSetSolution() first");
665   if (!ts->type_name) {
666     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
667   }
668   ierr = (*ts->setup)(ts);CHKERRQ(ierr);
669   ts->setupcalled = 1;
670   PetscFunctionReturn(0);
671 }
672 
673 #undef __FUNC__
674 #define __FUNC__ "TSDestroy"
675 /*@C
676    TSDestroy - Destroys the timestepper context that was created
677    with TSCreate().
678 
679    Collective on TS
680 
681    Input Parameter:
682 .  ts - the TS context obtained from TSCreate()
683 
684    Level: beginner
685 
686 .keywords: TS, timestepper, destroy
687 
688 .seealso: TSCreate(), TSSetUp(), TSSolve()
689 @*/
690 int TSDestroy(TS ts)
691 {
692   int ierr;
693 
694   PetscFunctionBegin;
695   PetscValidHeaderSpecific(ts,TS_COOKIE);
696   if (--ts->refct > 0) PetscFunctionReturn(0);
697 
698   /* if memory was published with AMS then destroy it */
699   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
700 
701   if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);}
702   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
703   ierr = (*(ts)->destroy)(ts);CHKERRQ(ierr);
704   PLogObjectDestroy((PetscObject)ts);
705   PetscHeaderDestroy((PetscObject)ts);
706   PetscFunctionReturn(0);
707 }
708 
709 #undef __FUNC__
710 #define __FUNC__ "TSGetSNES"
711 /*@C
712    TSGetSNES - Returns the SNES (nonlinear solver) associated with
713    a TS (timestepper) context. Valid only for nonlinear problems.
714 
715    Not Collective, but SNES is parallel if TS is parallel
716 
717    Input Parameter:
718 .  ts - the TS context obtained from TSCreate()
719 
720    Output Parameter:
721 .  snes - the nonlinear solver context
722 
723    Notes:
724    The user can then directly manipulate the SNES context to set various
725    options, etc.  Likewise, the user can then extract and manipulate the
726    SLES, KSP, and PC contexts as well.
727 
728    TSGetSNES() does not work for integrators that do not use SNES; in
729    this case TSGetSNES() returns PETSC_NULL in snes.
730 
731    Level: beginner
732 
733 .keywords: timestep, get, SNES
734 @*/
735 int TSGetSNES(TS ts,SNES *snes)
736 {
737   PetscFunctionBegin;
738   PetscValidHeaderSpecific(ts,TS_COOKIE);
739   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Nonlinear only; use TSGetSLES()");
740   *snes = ts->snes;
741   PetscFunctionReturn(0);
742 }
743 
744 #undef __FUNC__
745 #define __FUNC__ "TSGetSLES"
746 /*@C
747    TSGetSLES - Returns the SLES (linear solver) associated with
748    a TS (timestepper) context.
749 
750    Not Collective, but SLES is parallel if TS is parallel
751 
752    Input Parameter:
753 .  ts - the TS context obtained from TSCreate()
754 
755    Output Parameter:
756 .  sles - the nonlinear solver context
757 
758    Notes:
759    The user can then directly manipulate the SLES context to set various
760    options, etc.  Likewise, the user can then extract and manipulate the
761    KSP and PC contexts as well.
762 
763    TSGetSLES() does not work for integrators that do not use SLES;
764    in this case TSGetSLES() returns PETSC_NULL in sles.
765 
766    Level: beginner
767 
768 .keywords: timestep, get, SLES
769 @*/
770 int TSGetSLES(TS ts,SLES *sles)
771 {
772   PetscFunctionBegin;
773   PetscValidHeaderSpecific(ts,TS_COOKIE);
774   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Linear only; use TSGetSNES()");
775   *sles = ts->sles;
776   PetscFunctionReturn(0);
777 }
778 
779 /* ----------- Routines to set solver parameters ---------- */
780 
781 #undef __FUNC__
782 #define __FUNC__ "TSSetDuration"
783 /*@
784    TSSetDuration - Sets the maximum number of timesteps to use and
785    maximum time for iteration.
786 
787    Collective on TS
788 
789    Input Parameters:
790 +  ts - the TS context obtained from TSCreate()
791 .  maxsteps - maximum number of iterations to use
792 -  maxtime - final time to iterate to
793 
794    Options Database Keys:
795 .  -ts_max_steps <maxsteps> - Sets maxsteps
796 .  -ts_max_time <maxtime> - Sets maxtime
797 
798    Notes:
799    The default maximum number of iterations is 5000. Default time is 5.0
800 
801    Level: intermediate
802 
803 .keywords: TS, timestep, set, maximum, iterations
804 @*/
805 int TSSetDuration(TS ts,int maxsteps,double maxtime)
806 {
807   PetscFunctionBegin;
808   PetscValidHeaderSpecific(ts,TS_COOKIE);
809   ts->max_steps = maxsteps;
810   ts->max_time  = maxtime;
811   PetscFunctionReturn(0);
812 }
813 
814 #undef __FUNC__
815 #define __FUNC__ "TSSetSolution"
816 /*@
817    TSSetSolution - Sets the initial solution vector
818    for use by the TS routines.
819 
820    Collective on TS and Vec
821 
822    Input Parameters:
823 +  ts - the TS context obtained from TSCreate()
824 -  x - the solution vector
825 
826    Level: beginner
827 
828 .keywords: TS, timestep, set, solution, initial conditions
829 @*/
830 int TSSetSolution(TS ts,Vec x)
831 {
832   PetscFunctionBegin;
833   PetscValidHeaderSpecific(ts,TS_COOKIE);
834   ts->vec_sol        = ts->vec_sol_always = x;
835   PetscFunctionReturn(0);
836 }
837 
838 /* ------------ Routines to set performance monitoring options ----------- */
839 
840 #undef __FUNC__
841 #define __FUNC__ "TSSetMonitor"
842 /*@C
843    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
844    timestep to display the iteration's  progress.
845 
846    Collective on TS
847 
848    Input Parameters:
849 +  ts - the TS context obtained from TSCreate()
850 .  func - monitoring routine
851 -  mctx - [optional] user-defined context for private data for the
852           monitor routine (may be PETSC_NULL)
853 
854    Calling sequence of func:
855 $    int func(TS ts,int steps,double time,Vec x,void *mctx)
856 
857 +    ts - the TS context
858 .    steps - iteration number
859 .    time - current timestep
860 .    x - current iterate
861 -    mctx - [optional] monitoring context
862 
863    Notes:
864    This routine adds an additional monitor to the list of monitors that
865    already has been loaded.
866 
867    Level: intermediate
868 
869 .keywords: TS, timestep, set, monitor
870 
871 .seealso: TSDefaultMonitor(), TSClearMonitor()
872 @*/
873 int TSSetMonitor(TS ts, int (*monitor)(TS,int,double,Vec,void*), void *mctx )
874 {
875   PetscFunctionBegin;
876   PetscValidHeaderSpecific(ts,TS_COOKIE);
877   if (ts->numbermonitors >= MAXTSMONITORS) {
878     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
879   }
880   ts->monitor[ts->numbermonitors]           = monitor;
881   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
882   PetscFunctionReturn(0);
883 }
884 
885 #undef __FUNC__
886 #define __FUNC__ "TSClearMonitor"
887 /*@C
888    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
889 
890    Collective on TS
891 
892    Input Parameters:
893 .  ts - the TS context obtained from TSCreate()
894 
895    Notes:
896    There is no way to remove a single, specific monitor.
897 
898    Level: intermediate
899 
900 .keywords: TS, timestep, set, monitor
901 
902 .seealso: TSDefaultMonitor(), TSSetMonitor()
903 @*/
904 int TSClearMonitor(TS ts)
905 {
906   PetscFunctionBegin;
907   PetscValidHeaderSpecific(ts,TS_COOKIE);
908   ts->numbermonitors = 0;
909   PetscFunctionReturn(0);
910 }
911 
912 #undef __FUNC__
913 #define __FUNC__ "TSDefaultMonitor"
914 int TSDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
915 {
916   int ierr;
917 
918   PetscFunctionBegin;
919   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,time);CHKERRQ(ierr);
920   PetscFunctionReturn(0);
921 }
922 
923 #undef __FUNC__
924 #define __FUNC__ "TSStep"
925 /*@
926    TSStep - Steps the requested number of timesteps.
927 
928    Collective on TS
929 
930    Input Parameter:
931 .  ts - the TS context obtained from TSCreate()
932 
933    Output Parameters:
934 +  steps - number of iterations until termination
935 -  time - time until termination
936 
937    Level: beginner
938 
939 .keywords: TS, timestep, solve
940 
941 .seealso: TSCreate(), TSSetUp(), TSDestroy()
942 @*/
943 int TSStep(TS ts,int *steps,double *time)
944 {
945   int ierr,flg;
946 
947   PetscFunctionBegin;
948   PetscValidHeaderSpecific(ts,TS_COOKIE);
949   if (!ts->setupcalled) {ierr = TSSetUp(ts);CHKERRQ(ierr);}
950   PLogEventBegin(TS_Step,ts,0,0,0);
951   ierr = (*(ts)->step)(ts,steps,time);CHKERRQ(ierr);
952   PLogEventEnd(TS_Step,ts,0,0,0);
953   ierr = OptionsHasName(ts->prefix,"-ts_view",&flg);CHKERRQ(ierr);
954   if (flg) {
955     ierr = TSView(ts,VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
956   }
957   PetscFunctionReturn(0);
958 }
959 
960 #undef __FUNC__
961 #define __FUNC__ "TSMonitor"
962 /*
963      Runs the user provided monitor routines, if they exists.
964 */
965 int TSMonitor(TS ts,int step,double time,Vec x)
966 {
967   int i,ierr,n = ts->numbermonitors;
968 
969   PetscFunctionBegin;
970   for ( i=0; i<n; i++ ) {
971     ierr = (*ts->monitor[i])(ts,step,time,x,ts->monitorcontext[i]);CHKERRQ(ierr);
972   }
973   PetscFunctionReturn(0);
974 }
975 
976 /* ------------------------------------------------------------------------*/
977 
978 #undef __FUNC__
979 #define __FUNC__ "TSLGMonitorCreate"
980 /*@C
981    TSLGMonitorCreate - Creates a line graph context for use with
982    TS to monitor convergence of preconditioned residual norms.
983 
984    Collective on TS
985 
986    Input Parameters:
987 +  host - the X display to open, or null for the local machine
988 .  label - the title to put in the title bar
989 .  x, y - the screen coordinates of the upper left coordinate of
990           the window
991 -  m, n - the screen width and height in pixels
992 
993    Output Parameter:
994 .  draw - the drawing context
995 
996    Options Database Key:
997 .  -ts_xmonitor - automatically sets line graph monitor
998 
999    Notes:
1000    Use TSLGMonitorDestroy() to destroy this line graph, not DrawLGDestroy().
1001 
1002    Level: intermediate
1003 
1004 .keywords: TS, monitor, line graph, residual, create
1005 
1006 .seealso: TSLGMonitorDestroy(), TSSetMonitor()
1007 @*/
1008 int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,
1009                        int n, DrawLG *draw)
1010 {
1011   Draw win;
1012   int  ierr;
1013 
1014   PetscFunctionBegin;
1015   ierr = DrawOpenX(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1016   ierr = DrawLGCreate(win,1,draw);CHKERRQ(ierr);
1017   ierr = DrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1018 
1019   PLogObjectParent(*draw,win);
1020   PetscFunctionReturn(0);
1021 }
1022 
1023 #undef __FUNC__
1024 #define __FUNC__ "TSLGMonitor"
1025 int TSLGMonitor(TS ts,int n,double time,Vec v,void *monctx)
1026 {
1027   DrawLG lg = (DrawLG) monctx;
1028   double x,y = time;
1029   int    ierr;
1030 
1031   PetscFunctionBegin;
1032   if (!n) {ierr = DrawLGReset(lg);CHKERRQ(ierr);}
1033   x = (double) n;
1034   ierr = DrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1035   if (n < 20 || (n % 5)) {
1036     ierr = DrawLGDraw(lg);CHKERRQ(ierr);
1037   }
1038   PetscFunctionReturn(0);
1039 }
1040 
1041 #undef __FUNC__
1042 #define __FUNC__ "TSLGMonitorDestroy"
1043 /*@C
1044    TSLGMonitorDestroy - Destroys a line graph context that was created
1045    with TSLGMonitorCreate().
1046 
1047    Collective on DrawLG
1048 
1049    Input Parameter:
1050 .  draw - the drawing context
1051 
1052    Level: intermediate
1053 
1054 .keywords: TS, monitor, line graph, destroy
1055 
1056 .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1057 @*/
1058 int TSLGMonitorDestroy(DrawLG drawlg)
1059 {
1060   Draw draw;
1061   int  ierr;
1062 
1063   PetscFunctionBegin;
1064   ierr = DrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1065   ierr = DrawDestroy(draw);CHKERRQ(ierr);
1066   ierr = DrawLGDestroy(drawlg);CHKERRQ(ierr);
1067   PetscFunctionReturn(0);
1068 }
1069 
1070 #undef __FUNC__
1071 #define __FUNC__ "TSGetTime"
1072 /*@
1073    TSGetTime - Gets the current time.
1074 
1075    Not Collective
1076 
1077    Input Parameter:
1078 .  ts - the TS context obtained from TSCreate()
1079 
1080    Output Parameter:
1081 .  t  - the current time
1082 
1083    Contributed by: Matthew Knepley
1084 
1085    Level: beginner
1086 
1087 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1088 
1089 .keywords: TS, get, time
1090 @*/
1091 int TSGetTime(TS ts, double* t)
1092 {
1093   PetscFunctionBegin;
1094   PetscValidHeaderSpecific(ts, TS_COOKIE);
1095   *t = ts->ptime;
1096   PetscFunctionReturn(0);
1097 }
1098 
1099 #undef __FUNC__
1100 #define __FUNC__ "TSGetProblemType"
1101 /*@C
1102    TSGetProblemType - Returns the problem type of a TS (timestepper) context.
1103 
1104    Not Collective
1105 
1106    Input Parameter:
1107 .  ts   - The TS context obtained from TSCreate()
1108 
1109    Output Parameter:
1110 .  type - The problem type, TS_LINEAR or TS_NONLINEAR
1111 
1112    Level: intermediate
1113 
1114    Contributed by: Matthew Knepley
1115 
1116 .keywords: ts, get, type
1117 
1118 @*/
1119 int TSGetProblemType(TS ts, TSProblemType *type)
1120 {
1121   PetscFunctionBegin;
1122   PetscValidHeaderSpecific(ts, TS_COOKIE);
1123   *type = ts->problem_type;
1124   PetscFunctionReturn(0);
1125 }
1126 
1127 #undef __FUNC__
1128 #define __FUNC__ "TSSetOptionsPrefix"
1129 /*@C
1130    TSSetOptionsPrefix - Sets the prefix used for searching for all
1131    TS options in the database.
1132 
1133    Collective on TS
1134 
1135    Input Parameter:
1136 +  ts     - The TS context
1137 -  prefix - The prefix to prepend to all option names
1138 
1139    Notes:
1140    A hyphen (-) must NOT be given at the beginning of the prefix name.
1141    The first character of all runtime options is AUTOMATICALLY the
1142    hyphen.
1143 
1144    Contributed by: Matthew Knepley
1145 
1146    Level: advanced
1147 
1148 .keywords: TS, set, options, prefix, database
1149 
1150 .seealso: TSSetFromOptions()
1151 
1152 @*/
1153 int TSSetOptionsPrefix(TS ts, char *prefix)
1154 {
1155   int ierr;
1156 
1157   PetscFunctionBegin;
1158   PetscValidHeaderSpecific(ts, TS_COOKIE);
1159   ierr = PetscObjectSetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1160   switch(ts->problem_type) {
1161     case TS_NONLINEAR:
1162       ierr = SNESSetOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr);
1163       break;
1164     case TS_LINEAR:
1165       ierr = SLESSetOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr);
1166       break;
1167   }
1168   PetscFunctionReturn(0);
1169 }
1170 
1171 
1172 #undef __FUNC__
1173 #define __FUNC__ "TSAppendOptionsPrefix"
1174 /*@C
1175    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1176    TS options in the database.
1177 
1178    Collective on TS
1179 
1180    Input Parameter:
1181 +  ts     - The TS context
1182 -  prefix - The prefix to prepend to all option names
1183 
1184    Notes:
1185    A hyphen (-) must NOT be given at the beginning of the prefix name.
1186    The first character of all runtime options is AUTOMATICALLY the
1187    hyphen.
1188 
1189    Contributed by: Matthew Knepley
1190 
1191    Level: advanced
1192 
1193 .keywords: TS, append, options, prefix, database
1194 
1195 .seealso: TSGetOptionsPrefix()
1196 
1197 @*/
1198 int TSAppendOptionsPrefix(TS ts, char *prefix)
1199 {
1200   int ierr;
1201 
1202   PetscFunctionBegin;
1203   PetscValidHeaderSpecific(ts, TS_COOKIE);
1204   ierr = PetscObjectAppendOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1205   switch(ts->problem_type) {
1206     case TS_NONLINEAR:
1207       ierr = SNESAppendOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr);
1208       break;
1209     case TS_LINEAR:
1210       ierr = SLESAppendOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr);
1211       break;
1212   }
1213   PetscFunctionReturn(0);
1214 }
1215 
1216 #undef __FUNC__
1217 #define __FUNC__ "TSGetOptionsPrefix"
1218 /*@C
1219    TSGetOptionsPrefix - Sets the prefix used for searching for all
1220    TS options in the database.
1221 
1222    Not Collective
1223 
1224    Input Parameter:
1225 .  ts - The TS context
1226 
1227    Output Parameter:
1228 .  prefix - A pointer to the prefix string used
1229 
1230    Contributed by: Matthew Knepley
1231 
1232    Notes: On the fortran side, the user should pass in a string 'prifix' of
1233    sufficient length to hold the prefix.
1234 
1235    Level: intermediate
1236 
1237 .keywords: TS, get, options, prefix, database
1238 
1239 .seealso: TSAppendOptionsPrefix()
1240 @*/
1241 int TSGetOptionsPrefix(TS ts, char **prefix)
1242 {
1243   int ierr;
1244 
1245   PetscFunctionBegin;
1246   PetscValidHeaderSpecific(ts, TS_COOKIE);
1247   ierr = PetscObjectGetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1248   PetscFunctionReturn(0);
1249 }
1250 
1251 #undef __FUNC__
1252 #define __FUNC__ "TSGetRHSMatrix"
1253 /*@C
1254    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1255 
1256    Not Collective, but parallel objects are returned if TS is parallel
1257 
1258    Input Parameter:
1259 .  ts  - The TS context obtained from TSCreate()
1260 
1261    Output Parameters:
1262 +  A   - The matrix A, where U_t = A(t) U
1263 .  M   - The preconditioner matrix, usually the same as A
1264 -  ctx - User-defined context for matrix evaluation routine
1265 
1266    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1267 
1268    Contributed by: Matthew Knepley
1269 
1270    Level: intermediate
1271 
1272 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1273 
1274 .keywords: TS, timestep, get, matrix
1275 
1276 @*/
1277 int TSGetRHSMatrix(TS ts, Mat *A, Mat *M, void **ctx)
1278 {
1279   PetscFunctionBegin;
1280   PetscValidHeaderSpecific(ts, TS_COOKIE);
1281   if (A)   *A = ts->A;
1282   if (M)   *M = ts->B;
1283   if (ctx) *ctx = ts->jacP;
1284   PetscFunctionReturn(0);
1285 }
1286 
1287 #undef __FUNC__
1288 #define __FUNC__ "TSGetRHSJacobian"
1289 /*@C
1290    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1291 
1292    Not Collective, but parallel objects are returned if TS is parallel
1293 
1294    Input Parameter:
1295 .  ts  - The TS context obtained from TSCreate()
1296 
1297    Output Parameters:
1298 +  J   - The Jacobian J of F, where U_t = F(U,t)
1299 .  M   - The preconditioner matrix, usually the same as J
1300 - ctx - User-defined context for Jacobian evaluation routine
1301 
1302    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1303 
1304    Contributed by: Matthew Knepley
1305 
1306    Level: intermediate
1307 
1308 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1309 
1310 .keywords: TS, timestep, get, matrix, Jacobian
1311 @*/
1312 int TSGetRHSJacobian(TS ts, Mat *J, Mat *M, void **ctx)
1313 {
1314   int ierr;
1315 
1316   PetscFunctionBegin;
1317   ierr = TSGetRHSMatrix(ts, J, M, ctx);CHKERRQ(ierr);
1318   PetscFunctionReturn(0);
1319 }
1320 
1321 /*MC
1322    TSRegister - Adds a method to the timestepping solver package.
1323 
1324    Synopsis:
1325 
1326    TSRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(TS))
1327 
1328    Not collective
1329 
1330    Input Parameters:
1331 +  name_solver - name of a new user-defined solver
1332 .  path - path (either absolute or relative) the library containing this solver
1333 .  name_create - name of routine to create method context
1334 -  routine_create - routine to create method context
1335 
1336    Notes:
1337    TSRegister() may be called multiple times to add several user-defined solvers.
1338 
1339    If dynamic libraries are used, then the fourth input argument (routine_create)
1340    is ignored.
1341 
1342    Sample usage:
1343 .vb
1344    TSRegister("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
1345               "MySolverCreate",MySolverCreate);
1346 .ve
1347 
1348    Then, your solver can be chosen with the procedural interface via
1349 $     TSSetType(ts,"my_solver")
1350    or at runtime via the option
1351 $     -ts_type my_solver
1352 
1353    Level: advanced
1354 
1355    $PETSC_ARCH, $PETSC_DIR, $PETSC_LDIR, and $BOPT occuring in pathname will be replaced with appropriate values.
1356 
1357 .keywords: TS, register
1358 
1359 .seealso: TSRegisterAll(), TSRegisterDestroy()
1360 M*/
1361 
1362 #undef __FUNC__
1363 #define __FUNC__ "TSRegister_Private"
1364 int TSRegister_Private(char *sname,char *path,char *name,int (*function)(TS))
1365 {
1366   char fullname[256];
1367   int  ierr;
1368 
1369   PetscFunctionBegin;
1370   ierr = FListConcat_Private(path,name,fullname); CHKERRQ(ierr);
1371   FListAdd_Private(&TSList,sname,fullname,        (int (*)(void*))function);
1372   PetscFunctionReturn(0);
1373 }
1374