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