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