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