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