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