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