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