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