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