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