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