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