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