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