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