xref: /petsc/src/ts/interface/ts.c (revision 8f1adeca01e0cf022e9bfde6d3952ec0985a8331)
1 
2 #include <petsc-private/tsimpl.h>        /*I "petscts.h"  I*/
3 #include <petscdmshell.h>
4 
5 /* Logging support */
6 PetscClassId  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   PetscBool      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_final_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  TSSetFromOptions(TS ts)
73 {
74   PetscBool      opt,flg;
75   PetscErrorCode ierr;
76   PetscViewer    monviewer;
77   char           monfilename[PETSC_MAX_PATH_LEN];
78   SNES           snes;
79   TSAdapt        adapt;
80 
81   PetscFunctionBegin;
82   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
83   ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr);
84     /* Handle TS type options */
85     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
86 
87     /* Handle generic TS options */
88     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
89     ierr = PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
90     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);CHKERRQ(ierr);
91     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&ts->time_step,PETSC_NULL);CHKERRQ(ierr);
92     opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time;
93     ierr = PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);CHKERRQ(ierr);
94     if (flg) {ierr = TSSetExactFinalTime(ts,opt);CHKERRQ(ierr);}
95     ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr);
96     ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr);
97     ierr = PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr);
98     ierr = PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,PETSC_NULL);CHKERRQ(ierr);
99     ierr = PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,PETSC_NULL);CHKERRQ(ierr);
100 
101     /* Monitor options */
102     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
103     if (flg) {
104       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr);
105       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
106     }
107     ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
108     if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);}
109 
110     opt  = PETSC_FALSE;
111     ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
112     if (opt) {
113       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
114     }
115     opt  = PETSC_FALSE;
116     ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
117     if (opt) {
118       void *ctx;
119       ierr = TSMonitorSolutionCreate(ts,PETSC_NULL,&ctx);CHKERRQ(ierr);
120       ierr = TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);CHKERRQ(ierr);
121     }
122     opt  = PETSC_FALSE;
123     ierr = PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
124     if (flg) {
125       PetscViewer ctx;
126       if (monfilename[0]) {
127         ierr = PetscViewerBinaryOpen(((PetscObject)ts)->comm,monfilename,FILE_MODE_WRITE,&ctx);CHKERRQ(ierr);
128       } else {
129         ctx = PETSC_VIEWER_BINARY_(((PetscObject)ts)->comm);
130       }
131       ierr = TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
132     }
133     opt  = PETSC_FALSE;
134     ierr = PetscOptionsString("-ts_monitor_solution_vtk","Save each time step to a binary file, use filename-%%03D.vts","TSMonitorSolutionVTK",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
135     if (flg) {
136       const char *ptr,*ptr2;
137       char *filetemplate;
138       if (!monfilename[0]) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
139       /* Do some cursory validation of the input. */
140       ierr = PetscStrstr(monfilename,"%",(char**)&ptr);CHKERRQ(ierr);
141       if (!ptr) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
142       for (ptr++ ; ptr && *ptr; ptr++) {
143         ierr = PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);CHKERRQ(ierr);
144         if (!ptr2 && (*ptr < '0' || '9' < *ptr)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03D.vts");
145         if (ptr2) break;
146       }
147       ierr = PetscStrallocpy(monfilename,&filetemplate);CHKERRQ(ierr);
148       ierr = TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);CHKERRQ(ierr);
149     }
150 
151     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
152     ierr = TSAdaptSetFromOptions(adapt);CHKERRQ(ierr);
153 
154     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
155     if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
156 
157     /* Handle specific TS options */
158     if (ts->ops->setfromoptions) {
159       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
160     }
161 
162     /* process any options handlers added with PetscObjectAddOptionsHandler() */
163     ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
164   ierr = PetscOptionsEnd();CHKERRQ(ierr);
165   PetscFunctionReturn(0);
166 }
167 
168 #undef __FUNCT__
169 #undef __FUNCT__
170 #define __FUNCT__ "TSComputeRHSJacobian"
171 /*@
172    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
173       set with TSSetRHSJacobian().
174 
175    Collective on TS and Vec
176 
177    Input Parameters:
178 +  ts - the TS context
179 .  t - current timestep
180 -  x - input vector
181 
182    Output Parameters:
183 +  A - Jacobian matrix
184 .  B - optional preconditioning matrix
185 -  flag - flag indicating matrix structure
186 
187    Notes:
188    Most users should not need to explicitly call this routine, as it
189    is used internally within the nonlinear solvers.
190 
191    See KSPSetOperators() for important information about setting the
192    flag parameter.
193 
194    Level: developer
195 
196 .keywords: SNES, compute, Jacobian, matrix
197 
198 .seealso:  TSSetRHSJacobian(), KSPSetOperators()
199 @*/
200 PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
201 {
202   PetscErrorCode ierr;
203   PetscInt       Xstate;
204   DM             dm;
205   TSDM           tsdm;
206   TSRHSJacobian  rhsjacobianfunc;
207   void           *ctx;
208   TSIJacobian    ijacobianfunc;
209 
210   PetscFunctionBegin;
211   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
212   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
213   PetscCheckSameComm(ts,1,X,3);
214   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
215   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
216   ierr = DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);CHKERRQ(ierr);
217   ierr = DMTSGetIJacobian(dm,&ijacobianfunc,PETSC_NULL);CHKERRQ(ierr);
218   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
219   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) {
220     *flg = ts->rhsjacobian.mstructure;
221     PetscFunctionReturn(0);
222   }
223 
224   if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
225 
226   if (rhsjacobianfunc) {
227     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
228     *flg = DIFFERENT_NONZERO_PATTERN;
229     PetscStackPush("TS user Jacobian function");
230     ierr = (*rhsjacobianfunc)(ts,t,X,A,B,flg,ctx);CHKERRQ(ierr);
231     PetscStackPop;
232     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
233     /* make sure user returned a correct Jacobian and preconditioner */
234     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
235     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
236   } else {
237     ierr = MatZeroEntries(*A);CHKERRQ(ierr);
238     if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);}
239     *flg = SAME_NONZERO_PATTERN;
240   }
241   ts->rhsjacobian.time = t;
242   ts->rhsjacobian.X = X;
243   ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr);
244   ts->rhsjacobian.mstructure = *flg;
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "TSComputeRHSFunction"
250 /*@
251    TSComputeRHSFunction - Evaluates the right-hand-side function.
252 
253    Collective on TS and Vec
254 
255    Input Parameters:
256 +  ts - the TS context
257 .  t - current time
258 -  x - state vector
259 
260    Output Parameter:
261 .  y - right hand side
262 
263    Note:
264    Most users should not need to explicitly call this routine, as it
265    is used internally within the nonlinear solvers.
266 
267    Level: developer
268 
269 .keywords: TS, compute
270 
271 .seealso: TSSetRHSFunction(), TSComputeIFunction()
272 @*/
273 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
274 {
275   PetscErrorCode ierr;
276 
277   TSRHSFunction rhsfunction;
278   TSIFunction   ifunction;
279   void          *ctx;
280   DM            dm;
281 
282   PetscFunctionBegin;
283   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
284   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
285   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
286   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
287   ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctx);CHKERRQ(ierr);
288   ierr = DMTSGetIFunction(dm,&ifunction,PETSC_NULL);CHKERRQ(ierr);
289 
290   if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
291 
292   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
293   if (rhsfunction) {
294     PetscStackPush("TS user right-hand-side function");
295     ierr = (*rhsfunction)(ts,t,x,y,ctx);CHKERRQ(ierr);
296     PetscStackPop;
297   } else {
298     ierr = VecZeroEntries(y);CHKERRQ(ierr);
299   }
300 
301   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
302   PetscFunctionReturn(0);
303 }
304 
305 #undef __FUNCT__
306 #define __FUNCT__ "TSGetRHSVec_Private"
307 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
308 {
309   Vec            F;
310   PetscErrorCode ierr;
311 
312   PetscFunctionBegin;
313   *Frhs = PETSC_NULL;
314   ierr = TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
315   if (!ts->Frhs) {
316     ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr);
317   }
318   *Frhs = ts->Frhs;
319   PetscFunctionReturn(0);
320 }
321 
322 #undef __FUNCT__
323 #define __FUNCT__ "TSGetRHSMats_Private"
324 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
325 {
326   Mat            A,B;
327   PetscErrorCode ierr;
328 
329   PetscFunctionBegin;
330   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
331   if (Arhs) {
332     if (!ts->Arhs) {
333       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
334     }
335     *Arhs = ts->Arhs;
336   }
337   if (Brhs) {
338     if (!ts->Brhs) {
339       ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
340     }
341     *Brhs = ts->Brhs;
342   }
343   PetscFunctionReturn(0);
344 }
345 
346 #undef __FUNCT__
347 #define __FUNCT__ "TSComputeIFunction"
348 /*@
349    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
350 
351    Collective on TS and Vec
352 
353    Input Parameters:
354 +  ts - the TS context
355 .  t - current time
356 .  X - state vector
357 .  Xdot - time derivative of state vector
358 -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
359 
360    Output Parameter:
361 .  Y - right hand side
362 
363    Note:
364    Most users should not need to explicitly call this routine, as it
365    is used internally within the nonlinear solvers.
366 
367    If the user did did not write their equations in implicit form, this
368    function recasts them in implicit form.
369 
370    Level: developer
371 
372 .keywords: TS, compute
373 
374 .seealso: TSSetIFunction(), TSComputeRHSFunction()
375 @*/
376 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
377 {
378   PetscErrorCode ierr;
379   TSIFunction    ifunction;
380   TSRHSFunction  rhsfunction;
381   void           *ctx;
382   DM             dm;
383 
384   PetscFunctionBegin;
385   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
386   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
387   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
388   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
389 
390   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
391   ierr = DMTSGetIFunction(dm,&ifunction,&ctx);CHKERRQ(ierr);
392   ierr = DMTSGetRHSFunction(dm,&rhsfunction,PETSC_NULL);CHKERRQ(ierr);
393 
394   if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
395 
396   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
397   if (ifunction) {
398     PetscStackPush("TS user implicit function");
399     ierr = (*ifunction)(ts,t,X,Xdot,Y,ctx);CHKERRQ(ierr);
400     PetscStackPop;
401   }
402   if (imex) {
403     if (!ifunction) {
404       ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);
405     }
406   } else if (rhsfunction) {
407     if (ifunction) {
408       Vec Frhs;
409       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
410       ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr);
411       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
412     } else {
413       ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr);
414       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
415     }
416   }
417   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
418   PetscFunctionReturn(0);
419 }
420 
421 #undef __FUNCT__
422 #define __FUNCT__ "TSComputeIJacobian"
423 /*@
424    TSComputeIJacobian - Evaluates the Jacobian of the DAE
425 
426    Collective on TS and Vec
427 
428    Input
429       Input Parameters:
430 +  ts - the TS context
431 .  t - current timestep
432 .  X - state vector
433 .  Xdot - time derivative of state vector
434 .  shift - shift to apply, see note below
435 -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
436 
437    Output Parameters:
438 +  A - Jacobian matrix
439 .  B - optional preconditioning matrix
440 -  flag - flag indicating matrix structure
441 
442    Notes:
443    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
444 
445    dF/dX + shift*dF/dXdot
446 
447    Most users should not need to explicitly call this routine, as it
448    is used internally within the nonlinear solvers.
449 
450    Level: developer
451 
452 .keywords: TS, compute, Jacobian, matrix
453 
454 .seealso:  TSSetIJacobian()
455 @*/
456 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
457 {
458   PetscInt Xstate, Xdotstate;
459   PetscErrorCode ierr;
460   TSIJacobian    ijacobian;
461   TSRHSJacobian  rhsjacobian;
462   DM             dm;
463   void           *ctx;
464 
465   PetscFunctionBegin;
466 
467   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
468   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
469   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
470   PetscValidPointer(A,6);
471   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
472   PetscValidPointer(B,7);
473   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
474   PetscValidPointer(flg,8);
475 
476   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
477   ierr = DMTSGetIJacobian(dm,&ijacobian,&ctx);CHKERRQ(ierr);
478   ierr = DMTSGetRHSJacobian(dm,&rhsjacobian,PETSC_NULL);CHKERRQ(ierr);
479 
480   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
481   ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr);
482   if (ts->ijacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->ijacobian.X == X && ts->ijacobian.Xstate == Xstate && ts->ijacobian.Xdot == Xdot && ts->ijacobian.Xdotstate == Xdotstate && ts->ijacobian.imex == imex))) {
483     *flg = ts->ijacobian.mstructure;
484     ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr);
485     PetscFunctionReturn(0);
486   }
487 
488   if (!rhsjacobian && !ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
489 
490   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
491   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
492   if (ijacobian) {
493     *flg = DIFFERENT_NONZERO_PATTERN;
494     PetscStackPush("TS user implicit Jacobian");
495     ierr = (*ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ctx);CHKERRQ(ierr);
496     PetscStackPop;
497     /* make sure user returned a correct Jacobian and preconditioner */
498     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
499     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
500   }
501   if (imex) {
502     if (!ijacobian) {  /* system was written as Xdot = F(t,X) */
503       ierr = MatZeroEntries(*A);CHKERRQ(ierr);
504       ierr = MatShift(*A,shift);CHKERRQ(ierr);
505       if (*A != *B) {
506         ierr = MatZeroEntries(*B);CHKERRQ(ierr);
507         ierr = MatShift(*B,shift);CHKERRQ(ierr);
508       }
509       *flg = SAME_PRECONDITIONER;
510     }
511   } else {
512     if (!ijacobian) {
513       ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr);
514       ierr = MatScale(*A,-1);CHKERRQ(ierr);
515       ierr = MatShift(*A,shift);CHKERRQ(ierr);
516       if (*A != *B) {
517         ierr = MatScale(*B,-1);CHKERRQ(ierr);
518         ierr = MatShift(*B,shift);CHKERRQ(ierr);
519       }
520     } else if (rhsjacobian) {
521       Mat Arhs,Brhs;
522       MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
523       ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
524       ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
525       axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
526       ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr);
527       if (*A != *B) {
528         ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr);
529       }
530       *flg = PetscMin(*flg,flg2);
531     }
532   }
533 
534   ts->ijacobian.time = t;
535   ts->ijacobian.X = X;
536   ts->ijacobian.Xdot = Xdot;
537   ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr);
538   ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr);
539   ts->ijacobian.shift = shift;
540   ts->ijacobian.imex = imex;
541   ts->ijacobian.mstructure = *flg;
542   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
543   PetscFunctionReturn(0);
544 }
545 
546 #undef __FUNCT__
547 #define __FUNCT__ "TSSetRHSFunction"
548 /*@C
549     TSSetRHSFunction - Sets the routine for evaluating the function,
550     F(t,u), where U_t = F(t,u).
551 
552     Logically Collective on TS
553 
554     Input Parameters:
555 +   ts - the TS context obtained from TSCreate()
556 .   r - vector to put the computed right hand side (or PETSC_NULL to have it created)
557 .   f - routine for evaluating the right-hand-side function
558 -   ctx - [optional] user-defined context for private data for the
559           function evaluation routine (may be PETSC_NULL)
560 
561     Calling sequence of func:
562 $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
563 
564 +   t - current timestep
565 .   u - input vector
566 .   F - function vector
567 -   ctx - [optional] user-defined function context
568 
569     Level: beginner
570 
571 .keywords: TS, timestep, set, right-hand-side, function
572 
573 .seealso: TSSetRHSJacobian(), TSSetIJacobian()
574 @*/
575 PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
576 {
577   PetscErrorCode ierr;
578   SNES           snes;
579   Vec            ralloc = PETSC_NULL;
580   DM             dm;
581 
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
584   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
585 
586   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
587   ierr = DMTSSetRHSFunction(dm,f,ctx);CHKERRQ(ierr);
588   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
589   if (!r && !ts->dm && ts->vec_sol) {
590     ierr = VecDuplicate(ts->vec_sol,&ralloc);CHKERRQ(ierr);
591     r = ralloc;
592   }
593   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
594   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "TSSetRHSJacobian"
600 /*@C
601    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
602    where U_t = F(U,t), as well as the location to store the matrix.
603 
604    Logically Collective on TS
605 
606    Input Parameters:
607 +  ts  - the TS context obtained from TSCreate()
608 .  A   - Jacobian matrix
609 .  B   - preconditioner matrix (usually same as A)
610 .  f   - the Jacobian evaluation routine
611 -  ctx - [optional] user-defined context for private data for the
612          Jacobian evaluation routine (may be PETSC_NULL)
613 
614    Calling sequence of func:
615 $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
616 
617 +  t - current timestep
618 .  u - input vector
619 .  A - matrix A, where U_t = A(t)u
620 .  B - preconditioner matrix, usually the same as A
621 .  flag - flag indicating information about the preconditioner matrix
622           structure (same as flag in KSPSetOperators())
623 -  ctx - [optional] user-defined context for matrix evaluation routine
624 
625    Notes:
626    See KSPSetOperators() for important information about setting the flag
627    output parameter in the routine func().  Be sure to read this information!
628 
629    The routine func() takes Mat * as the matrix arguments rather than Mat.
630    This allows the matrix evaluation routine to replace A and/or B with a
631    completely new matrix structure (not just different matrix elements)
632    when appropriate, for instance, if the nonzero structure is changing
633    throughout the global iterations.
634 
635    Level: beginner
636 
637 .keywords: TS, timestep, set, right-hand-side, Jacobian
638 
639 .seealso: SNESDefaultComputeJacobianColor(), TSSetRHSFunction()
640 
641 @*/
642 PetscErrorCode  TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
643 {
644   PetscErrorCode ierr;
645   SNES           snes;
646   DM             dm;
647   TSIJacobian    ijacobian;
648 
649   PetscFunctionBegin;
650   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
651   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
652   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
653   if (A) PetscCheckSameComm(ts,1,A,2);
654   if (B) PetscCheckSameComm(ts,1,B,3);
655 
656   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
657   ierr = DMTSSetRHSJacobian(dm,f,ctx);CHKERRQ(ierr);
658   ierr = DMTSGetIJacobian(dm,&ijacobian,PETSC_NULL);CHKERRQ(ierr);
659 
660   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
661   if (!ijacobian) {
662     ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
663   }
664   if (A) {
665     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
666     ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
667     ts->Arhs = A;
668   }
669   if (B) {
670     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
671     ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
672     ts->Brhs = B;
673   }
674   PetscFunctionReturn(0);
675 }
676 
677 
678 #undef __FUNCT__
679 #define __FUNCT__ "TSSetIFunction"
680 /*@C
681    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
682 
683    Logically Collective on TS
684 
685    Input Parameters:
686 +  ts  - the TS context obtained from TSCreate()
687 .  r   - vector to hold the residual (or PETSC_NULL to have it created internally)
688 .  f   - the function evaluation routine
689 -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
690 
691    Calling sequence of f:
692 $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
693 
694 +  t   - time at step/stage being solved
695 .  u   - state vector
696 .  u_t - time derivative of state vector
697 .  F   - function vector
698 -  ctx - [optional] user-defined context for matrix evaluation routine
699 
700    Important:
701    The user MUST call either this routine, TSSetRHSFunction().  This routine must be used when not solving an ODE, for example a DAE.
702 
703    Level: beginner
704 
705 .keywords: TS, timestep, set, DAE, Jacobian
706 
707 .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
708 @*/
709 PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
710 {
711   PetscErrorCode ierr;
712   SNES           snes;
713   Vec            resalloc = PETSC_NULL;
714   DM             dm;
715 
716   PetscFunctionBegin;
717   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
718   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
719 
720   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
721   ierr = DMTSSetIFunction(dm,f,ctx);CHKERRQ(ierr);
722 
723   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
724   if (!res && !ts->dm && ts->vec_sol) {
725     ierr = VecDuplicate(ts->vec_sol,&resalloc);CHKERRQ(ierr);
726     res = resalloc;
727   }
728   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
729   ierr = VecDestroy(&resalloc);CHKERRQ(ierr);
730 
731   PetscFunctionReturn(0);
732 }
733 
734 #undef __FUNCT__
735 #define __FUNCT__ "TSGetIFunction"
736 /*@C
737    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
738 
739    Not Collective
740 
741    Input Parameter:
742 .  ts - the TS context
743 
744    Output Parameter:
745 +  r - vector to hold residual (or PETSC_NULL)
746 .  func - the function to compute residual (or PETSC_NULL)
747 -  ctx - the function context (or PETSC_NULL)
748 
749    Level: advanced
750 
751 .keywords: TS, nonlinear, get, function
752 
753 .seealso: TSSetIFunction(), SNESGetFunction()
754 @*/
755 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
756 {
757   PetscErrorCode ierr;
758   SNES snes;
759   DM   dm;
760 
761   PetscFunctionBegin;
762   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
763   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
764   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
765   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
766   ierr = DMTSGetIFunction(dm,func,ctx);CHKERRQ(ierr);
767   PetscFunctionReturn(0);
768 }
769 
770 #undef __FUNCT__
771 #define __FUNCT__ "TSGetRHSFunction"
772 /*@C
773    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
774 
775    Not Collective
776 
777    Input Parameter:
778 .  ts - the TS context
779 
780    Output Parameter:
781 +  r - vector to hold computed right hand side (or PETSC_NULL)
782 .  func - the function to compute right hand side (or PETSC_NULL)
783 -  ctx - the function context (or PETSC_NULL)
784 
785    Level: advanced
786 
787 .keywords: TS, nonlinear, get, function
788 
789 .seealso: TSSetRhsfunction(), SNESGetFunction()
790 @*/
791 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
792 {
793   PetscErrorCode ierr;
794   SNES snes;
795   DM   dm;
796 
797   PetscFunctionBegin;
798   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
799   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
800   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
801   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
802   ierr = DMTSGetRHSFunction(dm,func,ctx);CHKERRQ(ierr);
803   PetscFunctionReturn(0);
804 }
805 
806 #undef __FUNCT__
807 #define __FUNCT__ "TSSetIJacobian"
808 /*@C
809    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
810         you provided with TSSetIFunction().
811 
812    Logically Collective on TS
813 
814    Input Parameters:
815 +  ts  - the TS context obtained from TSCreate()
816 .  A   - Jacobian matrix
817 .  B   - preconditioning matrix for A (may be same as A)
818 .  f   - the Jacobian evaluation routine
819 -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
820 
821    Calling sequence of f:
822 $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
823 
824 +  t    - time at step/stage being solved
825 .  U    - state vector
826 .  U_t  - time derivative of state vector
827 .  a    - shift
828 .  A    - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
829 .  B    - preconditioning matrix for A, may be same as A
830 .  flag - flag indicating information about the preconditioner matrix
831           structure (same as flag in KSPSetOperators())
832 -  ctx  - [optional] user-defined context for matrix evaluation routine
833 
834    Notes:
835    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
836 
837    The matrix dF/dU + a*dF/dU_t you provide turns out to be
838    the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
839    The time integrator internally approximates U_t by W+a*U where the positive "shift"
840    a and vector W depend on the integration method, step size, and past states. For example with
841    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
842    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
843 
844    Level: beginner
845 
846 .keywords: TS, timestep, DAE, Jacobian
847 
848 .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESDefaultComputeJacobianColor(), SNESDefaultComputeJacobian()
849 
850 @*/
851 PetscErrorCode  TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
852 {
853   PetscErrorCode ierr;
854   SNES           snes;
855   DM             dm;
856 
857   PetscFunctionBegin;
858   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
859   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
860   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
861   if (A) PetscCheckSameComm(ts,1,A,2);
862   if (B) PetscCheckSameComm(ts,1,B,3);
863 
864   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
865   ierr = DMTSSetIJacobian(dm,f,ctx);CHKERRQ(ierr);
866 
867   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
868   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
869   PetscFunctionReturn(0);
870 }
871 
872 #undef __FUNCT__
873 #define __FUNCT__ "TSView"
874 /*@C
875     TSView - Prints the TS data structure.
876 
877     Collective on TS
878 
879     Input Parameters:
880 +   ts - the TS context obtained from TSCreate()
881 -   viewer - visualization context
882 
883     Options Database Key:
884 .   -ts_view - calls TSView() at end of TSStep()
885 
886     Notes:
887     The available visualization contexts include
888 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
889 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
890          output where only the first processor opens
891          the file.  All other processors send their
892          data to the first processor to print.
893 
894     The user can open an alternative visualization context with
895     PetscViewerASCIIOpen() - output to a specified file.
896 
897     Level: beginner
898 
899 .keywords: TS, timestep, view
900 
901 .seealso: PetscViewerASCIIOpen()
902 @*/
903 PetscErrorCode  TSView(TS ts,PetscViewer viewer)
904 {
905   PetscErrorCode ierr;
906   const TSType   type;
907   PetscBool      iascii,isstring,isundials;
908 
909   PetscFunctionBegin;
910   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
911   if (!viewer) {
912     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
913   }
914   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
915   PetscCheckSameComm(ts,1,viewer,2);
916 
917   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
918   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
919   if (iascii) {
920     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr);
921     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
922     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
923     if (ts->problem_type == TS_NONLINEAR) {
924       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);CHKERRQ(ierr);
925       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr);
926     }
927     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);CHKERRQ(ierr);
928     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
929     if (ts->ops->view) {
930       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
931       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
932       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
933     }
934   } else if (isstring) {
935     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
936     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
937   }
938   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
939   ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
940   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
941   PetscFunctionReturn(0);
942 }
943 
944 
945 #undef __FUNCT__
946 #define __FUNCT__ "TSSetApplicationContext"
947 /*@
948    TSSetApplicationContext - Sets an optional user-defined context for
949    the timesteppers.
950 
951    Logically Collective on TS
952 
953    Input Parameters:
954 +  ts - the TS context obtained from TSCreate()
955 -  usrP - optional user context
956 
957    Level: intermediate
958 
959 .keywords: TS, timestep, set, application, context
960 
961 .seealso: TSGetApplicationContext()
962 @*/
963 PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
964 {
965   PetscFunctionBegin;
966   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
967   ts->user = usrP;
968   PetscFunctionReturn(0);
969 }
970 
971 #undef __FUNCT__
972 #define __FUNCT__ "TSGetApplicationContext"
973 /*@
974     TSGetApplicationContext - Gets the user-defined context for the
975     timestepper.
976 
977     Not Collective
978 
979     Input Parameter:
980 .   ts - the TS context obtained from TSCreate()
981 
982     Output Parameter:
983 .   usrP - user context
984 
985     Level: intermediate
986 
987 .keywords: TS, timestep, get, application, context
988 
989 .seealso: TSSetApplicationContext()
990 @*/
991 PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
992 {
993   PetscFunctionBegin;
994   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
995   *(void**)usrP = ts->user;
996   PetscFunctionReturn(0);
997 }
998 
999 #undef __FUNCT__
1000 #define __FUNCT__ "TSGetTimeStepNumber"
1001 /*@
1002    TSGetTimeStepNumber - Gets the number of time steps completed.
1003 
1004    Not Collective
1005 
1006    Input Parameter:
1007 .  ts - the TS context obtained from TSCreate()
1008 
1009    Output Parameter:
1010 .  iter - number of steps completed so far
1011 
1012    Level: intermediate
1013 
1014 .keywords: TS, timestep, get, iteration, number
1015 .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStep()
1016 @*/
1017 PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
1018 {
1019   PetscFunctionBegin;
1020   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1021   PetscValidIntPointer(iter,2);
1022   *iter = ts->steps;
1023   PetscFunctionReturn(0);
1024 }
1025 
1026 #undef __FUNCT__
1027 #define __FUNCT__ "TSSetInitialTimeStep"
1028 /*@
1029    TSSetInitialTimeStep - Sets the initial timestep to be used,
1030    as well as the initial time.
1031 
1032    Logically Collective on TS
1033 
1034    Input Parameters:
1035 +  ts - the TS context obtained from TSCreate()
1036 .  initial_time - the initial time
1037 -  time_step - the size of the timestep
1038 
1039    Level: intermediate
1040 
1041 .seealso: TSSetTimeStep(), TSGetTimeStep()
1042 
1043 .keywords: TS, set, initial, timestep
1044 @*/
1045 PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
1046 {
1047   PetscErrorCode ierr;
1048 
1049   PetscFunctionBegin;
1050   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1051   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
1052   ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr);
1053   PetscFunctionReturn(0);
1054 }
1055 
1056 #undef __FUNCT__
1057 #define __FUNCT__ "TSSetTimeStep"
1058 /*@
1059    TSSetTimeStep - Allows one to reset the timestep at any time,
1060    useful for simple pseudo-timestepping codes.
1061 
1062    Logically Collective on TS
1063 
1064    Input Parameters:
1065 +  ts - the TS context obtained from TSCreate()
1066 -  time_step - the size of the timestep
1067 
1068    Level: intermediate
1069 
1070 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1071 
1072 .keywords: TS, set, timestep
1073 @*/
1074 PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
1075 {
1076   PetscFunctionBegin;
1077   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1078   PetscValidLogicalCollectiveReal(ts,time_step,2);
1079   ts->time_step = time_step;
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 #undef __FUNCT__
1084 #define __FUNCT__ "TSSetExactFinalTime"
1085 /*@
1086    TSSetExactFinalTime - Determines whether to interpolate solution to the
1087       exact final time requested by the user or just returns it at the final time
1088       it computed.
1089 
1090   Logically Collective on TS
1091 
1092    Input Parameter:
1093 +   ts - the time-step context
1094 -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
1095 
1096    Level: beginner
1097 
1098 .seealso: TSSetDuration()
1099 @*/
1100 PetscErrorCode  TSSetExactFinalTime(TS ts,PetscBool flg)
1101 {
1102 
1103   PetscFunctionBegin;
1104   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1105   PetscValidLogicalCollectiveBool(ts,flg,2);
1106   ts->exact_final_time = flg;
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 #undef __FUNCT__
1111 #define __FUNCT__ "TSGetTimeStep"
1112 /*@
1113    TSGetTimeStep - Gets the current timestep size.
1114 
1115    Not Collective
1116 
1117    Input Parameter:
1118 .  ts - the TS context obtained from TSCreate()
1119 
1120    Output Parameter:
1121 .  dt - the current timestep size
1122 
1123    Level: intermediate
1124 
1125 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1126 
1127 .keywords: TS, get, timestep
1128 @*/
1129 PetscErrorCode  TSGetTimeStep(TS ts,PetscReal* dt)
1130 {
1131   PetscFunctionBegin;
1132   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1133   PetscValidRealPointer(dt,2);
1134   *dt = ts->time_step;
1135   PetscFunctionReturn(0);
1136 }
1137 
1138 #undef __FUNCT__
1139 #define __FUNCT__ "TSGetSolution"
1140 /*@
1141    TSGetSolution - Returns the solution at the present timestep. It
1142    is valid to call this routine inside the function that you are evaluating
1143    in order to move to the new timestep. This vector not changed until
1144    the solution at the next timestep has been calculated.
1145 
1146    Not Collective, but Vec returned is parallel if TS is parallel
1147 
1148    Input Parameter:
1149 .  ts - the TS context obtained from TSCreate()
1150 
1151    Output Parameter:
1152 .  v - the vector containing the solution
1153 
1154    Level: intermediate
1155 
1156 .seealso: TSGetTimeStep()
1157 
1158 .keywords: TS, timestep, get, solution
1159 @*/
1160 PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1161 {
1162   PetscFunctionBegin;
1163   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1164   PetscValidPointer(v,2);
1165   *v = ts->vec_sol;
1166   PetscFunctionReturn(0);
1167 }
1168 
1169 /* ----- Routines to initialize and destroy a timestepper ---- */
1170 #undef __FUNCT__
1171 #define __FUNCT__ "TSSetProblemType"
1172 /*@
1173   TSSetProblemType - Sets the type of problem to be solved.
1174 
1175   Not collective
1176 
1177   Input Parameters:
1178 + ts   - The TS
1179 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1180 .vb
1181          U_t = A U
1182          U_t = A(t) U
1183          U_t = F(t,U)
1184 .ve
1185 
1186    Level: beginner
1187 
1188 .keywords: TS, problem type
1189 .seealso: TSSetUp(), TSProblemType, TS
1190 @*/
1191 PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1192 {
1193   PetscErrorCode ierr;
1194 
1195   PetscFunctionBegin;
1196   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1197   ts->problem_type = type;
1198   if (type == TS_LINEAR) {
1199     SNES snes;
1200     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1201     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
1202   }
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 #undef __FUNCT__
1207 #define __FUNCT__ "TSGetProblemType"
1208 /*@C
1209   TSGetProblemType - Gets the type of problem to be solved.
1210 
1211   Not collective
1212 
1213   Input Parameter:
1214 . ts   - The TS
1215 
1216   Output Parameter:
1217 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1218 .vb
1219          M U_t = A U
1220          M(t) U_t = A(t) U
1221          U_t = F(t,U)
1222 .ve
1223 
1224    Level: beginner
1225 
1226 .keywords: TS, problem type
1227 .seealso: TSSetUp(), TSProblemType, TS
1228 @*/
1229 PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1230 {
1231   PetscFunctionBegin;
1232   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1233   PetscValidIntPointer(type,2);
1234   *type = ts->problem_type;
1235   PetscFunctionReturn(0);
1236 }
1237 
1238 #undef __FUNCT__
1239 #define __FUNCT__ "TSSetUp"
1240 /*@
1241    TSSetUp - Sets up the internal data structures for the later use
1242    of a timestepper.
1243 
1244    Collective on TS
1245 
1246    Input Parameter:
1247 .  ts - the TS context obtained from TSCreate()
1248 
1249    Notes:
1250    For basic use of the TS solvers the user need not explicitly call
1251    TSSetUp(), since these actions will automatically occur during
1252    the call to TSStep().  However, if one wishes to control this
1253    phase separately, TSSetUp() should be called after TSCreate()
1254    and optional routines of the form TSSetXXX(), but before TSStep().
1255 
1256    Level: advanced
1257 
1258 .keywords: TS, timestep, setup
1259 
1260 .seealso: TSCreate(), TSStep(), TSDestroy()
1261 @*/
1262 PetscErrorCode  TSSetUp(TS ts)
1263 {
1264   PetscErrorCode ierr;
1265   DM             dm;
1266   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
1267   PetscErrorCode (*jac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
1268   TSIJacobian    ijac;
1269   TSRHSJacobian  rhsjac;
1270 
1271   PetscFunctionBegin;
1272   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1273   if (ts->setupcalled) PetscFunctionReturn(0);
1274 
1275   if (!((PetscObject)ts)->type_name) {
1276     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1277   }
1278   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1279 
1280   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1281 
1282   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
1283 
1284   if (ts->ops->setup) {
1285     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1286   }
1287 
1288   /* in the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
1289    to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
1290    */
1291   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1292   ierr = DMSNESGetFunction(dm,&func,PETSC_NULL);CHKERRQ(ierr);
1293   if (!func) {
1294     ierr =DMSNESSetFunction(dm,SNESTSFormFunction,ts);CHKERRQ(ierr);
1295   }
1296   /* if the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
1297      Otherwise, the SNES will use coloring internally to form the Jacobian.
1298    */
1299   ierr = DMSNESGetJacobian(dm,&jac,PETSC_NULL);CHKERRQ(ierr);
1300   ierr = DMTSGetIJacobian(dm,&ijac,PETSC_NULL);CHKERRQ(ierr);
1301   ierr = DMTSGetRHSJacobian(dm,&rhsjac,PETSC_NULL);CHKERRQ(ierr);
1302   if (!jac && (ijac || rhsjac)) {
1303     ierr = DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);CHKERRQ(ierr);
1304   }
1305   ts->setupcalled = PETSC_TRUE;
1306   PetscFunctionReturn(0);
1307 }
1308 
1309 #undef __FUNCT__
1310 #define __FUNCT__ "TSReset"
1311 /*@
1312    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1313 
1314    Collective on TS
1315 
1316    Input Parameter:
1317 .  ts - the TS context obtained from TSCreate()
1318 
1319    Level: beginner
1320 
1321 .keywords: TS, timestep, reset
1322 
1323 .seealso: TSCreate(), TSSetup(), TSDestroy()
1324 @*/
1325 PetscErrorCode  TSReset(TS ts)
1326 {
1327   PetscErrorCode ierr;
1328 
1329   PetscFunctionBegin;
1330   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1331   if (ts->ops->reset) {
1332     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1333   }
1334   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
1335   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
1336   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1337   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
1338   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1339   ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
1340   ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
1341   ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);
1342   ts->setupcalled = PETSC_FALSE;
1343   PetscFunctionReturn(0);
1344 }
1345 
1346 #undef __FUNCT__
1347 #define __FUNCT__ "TSDestroy"
1348 /*@
1349    TSDestroy - Destroys the timestepper context that was created
1350    with TSCreate().
1351 
1352    Collective on TS
1353 
1354    Input Parameter:
1355 .  ts - the TS context obtained from TSCreate()
1356 
1357    Level: beginner
1358 
1359 .keywords: TS, timestepper, destroy
1360 
1361 .seealso: TSCreate(), TSSetUp(), TSSolve()
1362 @*/
1363 PetscErrorCode  TSDestroy(TS *ts)
1364 {
1365   PetscErrorCode ierr;
1366 
1367   PetscFunctionBegin;
1368   if (!*ts) PetscFunctionReturn(0);
1369   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
1370   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1371 
1372   ierr = TSReset((*ts));CHKERRQ(ierr);
1373 
1374   /* if memory was published with AMS then destroy it */
1375   ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr);
1376   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
1377 
1378   ierr = TSAdaptDestroy(&(*ts)->adapt);CHKERRQ(ierr);
1379   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
1380   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
1381   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
1382 
1383   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1384   PetscFunctionReturn(0);
1385 }
1386 
1387 #undef __FUNCT__
1388 #define __FUNCT__ "TSGetSNES"
1389 /*@
1390    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1391    a TS (timestepper) context. Valid only for nonlinear problems.
1392 
1393    Not Collective, but SNES is parallel if TS is parallel
1394 
1395    Input Parameter:
1396 .  ts - the TS context obtained from TSCreate()
1397 
1398    Output Parameter:
1399 .  snes - the nonlinear solver context
1400 
1401    Notes:
1402    The user can then directly manipulate the SNES context to set various
1403    options, etc.  Likewise, the user can then extract and manipulate the
1404    KSP, KSP, and PC contexts as well.
1405 
1406    TSGetSNES() does not work for integrators that do not use SNES; in
1407    this case TSGetSNES() returns PETSC_NULL in snes.
1408 
1409    Level: beginner
1410 
1411 .keywords: timestep, get, SNES
1412 @*/
1413 PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1414 {
1415   PetscErrorCode ierr;
1416 
1417   PetscFunctionBegin;
1418   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1419   PetscValidPointer(snes,2);
1420   if (!ts->snes) {
1421     ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
1422     ierr = SNESSetFunction(ts->snes,PETSC_NULL,SNESTSFormFunction,ts);CHKERRQ(ierr);
1423     ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr);
1424     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
1425     if (ts->dm) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
1426     if (ts->problem_type == TS_LINEAR) {
1427       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
1428     }
1429   }
1430   *snes = ts->snes;
1431   PetscFunctionReturn(0);
1432 }
1433 
1434 #undef __FUNCT__
1435 #define __FUNCT__ "TSGetKSP"
1436 /*@
1437    TSGetKSP - Returns the KSP (linear solver) associated with
1438    a TS (timestepper) context.
1439 
1440    Not Collective, but KSP is parallel if TS is parallel
1441 
1442    Input Parameter:
1443 .  ts - the TS context obtained from TSCreate()
1444 
1445    Output Parameter:
1446 .  ksp - the nonlinear solver context
1447 
1448    Notes:
1449    The user can then directly manipulate the KSP context to set various
1450    options, etc.  Likewise, the user can then extract and manipulate the
1451    KSP and PC contexts as well.
1452 
1453    TSGetKSP() does not work for integrators that do not use KSP;
1454    in this case TSGetKSP() returns PETSC_NULL in ksp.
1455 
1456    Level: beginner
1457 
1458 .keywords: timestep, get, KSP
1459 @*/
1460 PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1461 {
1462   PetscErrorCode ierr;
1463   SNES           snes;
1464 
1465   PetscFunctionBegin;
1466   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1467   PetscValidPointer(ksp,2);
1468   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1469   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1470   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1471   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
1472   PetscFunctionReturn(0);
1473 }
1474 
1475 /* ----------- Routines to set solver parameters ---------- */
1476 
1477 #undef __FUNCT__
1478 #define __FUNCT__ "TSGetDuration"
1479 /*@
1480    TSGetDuration - Gets the maximum number of timesteps to use and
1481    maximum time for iteration.
1482 
1483    Not Collective
1484 
1485    Input Parameters:
1486 +  ts       - the TS context obtained from TSCreate()
1487 .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1488 -  maxtime  - final time to iterate to, or PETSC_NULL
1489 
1490    Level: intermediate
1491 
1492 .keywords: TS, timestep, get, maximum, iterations, time
1493 @*/
1494 PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1495 {
1496   PetscFunctionBegin;
1497   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1498   if (maxsteps) {
1499     PetscValidIntPointer(maxsteps,2);
1500     *maxsteps = ts->max_steps;
1501   }
1502   if (maxtime) {
1503     PetscValidScalarPointer(maxtime,3);
1504     *maxtime  = ts->max_time;
1505   }
1506   PetscFunctionReturn(0);
1507 }
1508 
1509 #undef __FUNCT__
1510 #define __FUNCT__ "TSSetDuration"
1511 /*@
1512    TSSetDuration - Sets the maximum number of timesteps to use and
1513    maximum time for iteration.
1514 
1515    Logically Collective on TS
1516 
1517    Input Parameters:
1518 +  ts - the TS context obtained from TSCreate()
1519 .  maxsteps - maximum number of iterations to use
1520 -  maxtime - final time to iterate to
1521 
1522    Options Database Keys:
1523 .  -ts_max_steps <maxsteps> - Sets maxsteps
1524 .  -ts_final_time <maxtime> - Sets maxtime
1525 
1526    Notes:
1527    The default maximum number of iterations is 5000. Default time is 5.0
1528 
1529    Level: intermediate
1530 
1531 .keywords: TS, timestep, set, maximum, iterations
1532 
1533 .seealso: TSSetExactFinalTime()
1534 @*/
1535 PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1536 {
1537   PetscFunctionBegin;
1538   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1539   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1540   PetscValidLogicalCollectiveReal(ts,maxtime,2);
1541   if (maxsteps >= 0) ts->max_steps = maxsteps;
1542   if (maxtime != PETSC_DEFAULT) ts->max_time  = maxtime;
1543   PetscFunctionReturn(0);
1544 }
1545 
1546 #undef __FUNCT__
1547 #define __FUNCT__ "TSSetSolution"
1548 /*@
1549    TSSetSolution - Sets the initial solution vector
1550    for use by the TS routines.
1551 
1552    Logically Collective on TS and Vec
1553 
1554    Input Parameters:
1555 +  ts - the TS context obtained from TSCreate()
1556 -  x - the solution vector
1557 
1558    Level: beginner
1559 
1560 .keywords: TS, timestep, set, solution, initial conditions
1561 @*/
1562 PetscErrorCode  TSSetSolution(TS ts,Vec x)
1563 {
1564   PetscErrorCode ierr;
1565   DM             dm;
1566 
1567   PetscFunctionBegin;
1568   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1569   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1570   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
1571   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1572   ts->vec_sol = x;
1573   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1574   ierr = DMShellSetGlobalVector(dm,x);CHKERRQ(ierr);
1575   PetscFunctionReturn(0);
1576 }
1577 
1578 #undef __FUNCT__
1579 #define __FUNCT__ "TSSetPreStep"
1580 /*@C
1581   TSSetPreStep - Sets the general-purpose function
1582   called once at the beginning of each time step.
1583 
1584   Logically Collective on TS
1585 
1586   Input Parameters:
1587 + ts   - The TS context obtained from TSCreate()
1588 - func - The function
1589 
1590   Calling sequence of func:
1591 . func (TS ts);
1592 
1593   Level: intermediate
1594 
1595   Note:
1596   If a step is rejected, TSStep() will call this routine again before each attempt.
1597   The last completed time step number can be queried using TSGetTimeStepNumber(), the
1598   size of the step being attempted can be obtained using TSGetTimeStep().
1599 
1600 .keywords: TS, timestep
1601 .seealso: TSSetPreStage(), TSSetPostStep(), TSStep()
1602 @*/
1603 PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1604 {
1605   PetscFunctionBegin;
1606   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1607   ts->ops->prestep = func;
1608   PetscFunctionReturn(0);
1609 }
1610 
1611 #undef __FUNCT__
1612 #define __FUNCT__ "TSPreStep"
1613 /*@
1614   TSPreStep - Runs the user-defined pre-step function.
1615 
1616   Collective on TS
1617 
1618   Input Parameters:
1619 . ts   - The TS context obtained from TSCreate()
1620 
1621   Notes:
1622   TSPreStep() is typically used within time stepping implementations,
1623   so most users would not generally call this routine themselves.
1624 
1625   Level: developer
1626 
1627 .keywords: TS, timestep
1628 .seealso: TSSetPreStep(), TSPreStage(), TSPostStep()
1629 @*/
1630 PetscErrorCode  TSPreStep(TS ts)
1631 {
1632   PetscErrorCode ierr;
1633 
1634   PetscFunctionBegin;
1635   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1636   if (ts->ops->prestep) {
1637     PetscStackPush("TS PreStep function");
1638     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1639     PetscStackPop;
1640   }
1641   PetscFunctionReturn(0);
1642 }
1643 
1644 #undef __FUNCT__
1645 #define __FUNCT__ "TSSetPreStage"
1646 /*@C
1647   TSSetPreStage - Sets the general-purpose function
1648   called once at the beginning of each stage.
1649 
1650   Logically Collective on TS
1651 
1652   Input Parameters:
1653 + ts   - The TS context obtained from TSCreate()
1654 - func - The function
1655 
1656   Calling sequence of func:
1657 . PetscErrorCode func(TS ts, PetscReal stagetime);
1658 
1659   Level: intermediate
1660 
1661   Note:
1662   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
1663   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
1664   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
1665 
1666 .keywords: TS, timestep
1667 .seealso: TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
1668 @*/
1669 PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
1670 {
1671   PetscFunctionBegin;
1672   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1673   ts->ops->prestage = func;
1674   PetscFunctionReturn(0);
1675 }
1676 
1677 #undef __FUNCT__
1678 #define __FUNCT__ "TSPreStage"
1679 /*@
1680   TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()
1681 
1682   Collective on TS
1683 
1684   Input Parameters:
1685 . ts   - The TS context obtained from TSCreate()
1686 
1687   Notes:
1688   TSPreStage() is typically used within time stepping implementations,
1689   most users would not generally call this routine themselves.
1690 
1691   Level: developer
1692 
1693 .keywords: TS, timestep
1694 .seealso: TSSetPreStep(), TSPreStep(), TSPostStep()
1695 @*/
1696 PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
1697 {
1698   PetscErrorCode ierr;
1699 
1700   PetscFunctionBegin;
1701   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1702   if (ts->ops->prestage) {
1703     PetscStackPush("TS PreStage function");
1704     ierr = (*ts->ops->prestage)(ts,stagetime);CHKERRQ(ierr);
1705     PetscStackPop;
1706   }
1707   PetscFunctionReturn(0);
1708 }
1709 
1710 #undef __FUNCT__
1711 #define __FUNCT__ "TSSetPostStep"
1712 /*@C
1713   TSSetPostStep - Sets the general-purpose function
1714   called once at the end of each time step.
1715 
1716   Logically Collective on TS
1717 
1718   Input Parameters:
1719 + ts   - The TS context obtained from TSCreate()
1720 - func - The function
1721 
1722   Calling sequence of func:
1723 $ func (TS ts);
1724 
1725   Level: intermediate
1726 
1727 .keywords: TS, timestep
1728 .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime()
1729 @*/
1730 PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1731 {
1732   PetscFunctionBegin;
1733   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1734   ts->ops->poststep = func;
1735   PetscFunctionReturn(0);
1736 }
1737 
1738 #undef __FUNCT__
1739 #define __FUNCT__ "TSPostStep"
1740 /*@
1741   TSPostStep - Runs the user-defined post-step function.
1742 
1743   Collective on TS
1744 
1745   Input Parameters:
1746 . ts   - The TS context obtained from TSCreate()
1747 
1748   Notes:
1749   TSPostStep() is typically used within time stepping implementations,
1750   so most users would not generally call this routine themselves.
1751 
1752   Level: developer
1753 
1754 .keywords: TS, timestep
1755 @*/
1756 PetscErrorCode  TSPostStep(TS ts)
1757 {
1758   PetscErrorCode ierr;
1759 
1760   PetscFunctionBegin;
1761   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1762   if (ts->ops->poststep) {
1763     PetscStackPush("TS PostStep function");
1764     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1765     PetscStackPop;
1766   }
1767   PetscFunctionReturn(0);
1768 }
1769 
1770 /* ------------ Routines to set performance monitoring options ----------- */
1771 
1772 #undef __FUNCT__
1773 #define __FUNCT__ "TSMonitorSet"
1774 /*@C
1775    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1776    timestep to display the iteration's  progress.
1777 
1778    Logically Collective on TS
1779 
1780    Input Parameters:
1781 +  ts - the TS context obtained from TSCreate()
1782 .  monitor - monitoring routine
1783 .  mctx - [optional] user-defined context for private data for the
1784              monitor routine (use PETSC_NULL if no context is desired)
1785 -  monitordestroy - [optional] routine that frees monitor context
1786           (may be PETSC_NULL)
1787 
1788    Calling sequence of monitor:
1789 $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1790 
1791 +    ts - the TS context
1792 .    steps - iteration number
1793 .    time - current time
1794 .    x - current iterate
1795 -    mctx - [optional] monitoring context
1796 
1797    Notes:
1798    This routine adds an additional monitor to the list of monitors that
1799    already has been loaded.
1800 
1801    Fortran notes: Only a single monitor function can be set for each TS object
1802 
1803    Level: intermediate
1804 
1805 .keywords: TS, timestep, set, monitor
1806 
1807 .seealso: TSMonitorDefault(), TSMonitorCancel()
1808 @*/
1809 PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1810 {
1811   PetscFunctionBegin;
1812   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1813   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1814   ts->monitor[ts->numbermonitors]           = monitor;
1815   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1816   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1817   PetscFunctionReturn(0);
1818 }
1819 
1820 #undef __FUNCT__
1821 #define __FUNCT__ "TSMonitorCancel"
1822 /*@C
1823    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1824 
1825    Logically Collective on TS
1826 
1827    Input Parameters:
1828 .  ts - the TS context obtained from TSCreate()
1829 
1830    Notes:
1831    There is no way to remove a single, specific monitor.
1832 
1833    Level: intermediate
1834 
1835 .keywords: TS, timestep, set, monitor
1836 
1837 .seealso: TSMonitorDefault(), TSMonitorSet()
1838 @*/
1839 PetscErrorCode  TSMonitorCancel(TS ts)
1840 {
1841   PetscErrorCode ierr;
1842   PetscInt       i;
1843 
1844   PetscFunctionBegin;
1845   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1846   for (i=0; i<ts->numbermonitors; i++) {
1847     if (ts->mdestroy[i]) {
1848       ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
1849     }
1850   }
1851   ts->numbermonitors = 0;
1852   PetscFunctionReturn(0);
1853 }
1854 
1855 #undef __FUNCT__
1856 #define __FUNCT__ "TSMonitorDefault"
1857 /*@
1858    TSMonitorDefault - Sets the Default monitor
1859 
1860    Level: intermediate
1861 
1862 .keywords: TS, set, monitor
1863 
1864 .seealso: TSMonitorDefault(), TSMonitorSet()
1865 @*/
1866 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1867 {
1868   PetscErrorCode ierr;
1869   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1870 
1871   PetscFunctionBegin;
1872   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1873   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr);
1874   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1875   PetscFunctionReturn(0);
1876 }
1877 
1878 #undef __FUNCT__
1879 #define __FUNCT__ "TSSetRetainStages"
1880 /*@
1881    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
1882 
1883    Logically Collective on TS
1884 
1885    Input Argument:
1886 .  ts - time stepping context
1887 
1888    Output Argument:
1889 .  flg - PETSC_TRUE or PETSC_FALSE
1890 
1891    Level: intermediate
1892 
1893 .keywords: TS, set
1894 
1895 .seealso: TSInterpolate(), TSSetPostStep()
1896 @*/
1897 PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
1898 {
1899 
1900   PetscFunctionBegin;
1901   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1902   ts->retain_stages = flg;
1903   PetscFunctionReturn(0);
1904 }
1905 
1906 #undef __FUNCT__
1907 #define __FUNCT__ "TSInterpolate"
1908 /*@
1909    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
1910 
1911    Collective on TS
1912 
1913    Input Argument:
1914 +  ts - time stepping context
1915 -  t - time to interpolate to
1916 
1917    Output Argument:
1918 .  X - state at given time
1919 
1920    Notes:
1921    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
1922 
1923    Level: intermediate
1924 
1925    Developer Notes:
1926    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
1927 
1928 .keywords: TS, set
1929 
1930 .seealso: TSSetRetainStages(), TSSetPostStep()
1931 @*/
1932 PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
1933 {
1934   PetscErrorCode ierr;
1935 
1936   PetscFunctionBegin;
1937   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1938   if (t < ts->ptime - ts->time_step_prev || t > ts->ptime) SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested time %G not in last time steps [%G,%G]",t,ts->ptime-ts->time_step_prev,ts->ptime);
1939   if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
1940   ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr);
1941   PetscFunctionReturn(0);
1942 }
1943 
1944 #undef __FUNCT__
1945 #define __FUNCT__ "TSStep"
1946 /*@
1947    TSStep - Steps one time step
1948 
1949    Collective on TS
1950 
1951    Input Parameter:
1952 .  ts - the TS context obtained from TSCreate()
1953 
1954    Level: intermediate
1955 
1956    Notes:
1957    The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
1958    be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.
1959 
1960 .keywords: TS, timestep, solve
1961 
1962 .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage()
1963 @*/
1964 PetscErrorCode  TSStep(TS ts)
1965 {
1966   PetscReal      ptime_prev;
1967   PetscErrorCode ierr;
1968 
1969   PetscFunctionBegin;
1970   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1971   ierr = TSSetUp(ts);CHKERRQ(ierr);
1972 
1973   ts->reason = TS_CONVERGED_ITERATING;
1974 
1975   ptime_prev = ts->ptime;
1976   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1977   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
1978   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1979   ts->time_step_prev = ts->ptime - ptime_prev;
1980 
1981   if (ts->reason < 0) {
1982     if (ts->errorifstepfailed) {
1983       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
1984         SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
1985       } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
1986     }
1987   } else if (!ts->reason) {
1988     if (ts->steps >= ts->max_steps)
1989       ts->reason = TS_CONVERGED_ITS;
1990     else if (ts->ptime >= ts->max_time)
1991       ts->reason = TS_CONVERGED_TIME;
1992   }
1993 
1994   PetscFunctionReturn(0);
1995 }
1996 
1997 #undef __FUNCT__
1998 #define __FUNCT__ "TSEvaluateStep"
1999 /*@
2000    TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
2001 
2002    Collective on TS
2003 
2004    Input Arguments:
2005 +  ts - time stepping context
2006 .  order - desired order of accuracy
2007 -  done - whether the step was evaluated at this order (pass PETSC_NULL to generate an error if not available)
2008 
2009    Output Arguments:
2010 .  X - state at the end of the current step
2011 
2012    Level: advanced
2013 
2014    Notes:
2015    This function cannot be called until all stages have been evaluated.
2016    It is normally called by adaptive controllers before a step has been accepted and may also be called by the user after TSStep() has returned.
2017 
2018 .seealso: TSStep(), TSAdapt
2019 @*/
2020 PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done)
2021 {
2022   PetscErrorCode ierr;
2023 
2024   PetscFunctionBegin;
2025   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2026   PetscValidType(ts,1);
2027   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
2028   if (!ts->ops->evaluatestep) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
2029   ierr = (*ts->ops->evaluatestep)(ts,order,X,done);CHKERRQ(ierr);
2030   PetscFunctionReturn(0);
2031 }
2032 
2033 #undef __FUNCT__
2034 #define __FUNCT__ "TSSolve"
2035 /*@
2036    TSSolve - Steps the requested number of timesteps.
2037 
2038    Collective on TS
2039 
2040    Input Parameter:
2041 +  ts - the TS context obtained from TSCreate()
2042 -  x - the solution vector
2043 
2044    Output Parameter:
2045 .  ftime - time of the state vector x upon completion
2046 
2047    Level: beginner
2048 
2049    Notes:
2050    The final time returned by this function may be different from the time of the internally
2051    held state accessible by TSGetSolution() and TSGetTime() because the method may have
2052    stepped over the final time.
2053 
2054 .keywords: TS, timestep, solve
2055 
2056 .seealso: TSCreate(), TSSetSolution(), TSStep()
2057 @*/
2058 PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime)
2059 {
2060   PetscBool      flg;
2061   char           filename[PETSC_MAX_PATH_LEN];
2062   PetscViewer    viewer;
2063   PetscErrorCode ierr;
2064 
2065   PetscFunctionBegin;
2066   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2067   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
2068   if (ts->exact_final_time) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
2069     if (!ts->vec_sol || x == ts->vec_sol) {
2070       Vec y;
2071       ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
2072       ierr = TSSetSolution(ts,y);CHKERRQ(ierr);
2073       ierr = VecDestroy(&y);CHKERRQ(ierr); /* grant ownership */
2074     }
2075     ierr = VecCopy(x,ts->vec_sol);CHKERRQ(ierr);
2076   } else {
2077     ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
2078   }
2079   ierr = TSSetUp(ts);CHKERRQ(ierr);
2080   /* reset time step and iteration counters */
2081   ts->steps = 0;
2082   ts->ksp_its = 0;
2083   ts->snes_its = 0;
2084   ts->num_snes_failures = 0;
2085   ts->reject = 0;
2086   ts->reason = TS_CONVERGED_ITERATING;
2087 
2088   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
2089     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
2090     ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
2091     if (ftime) *ftime = ts->ptime;
2092   } else {
2093     /* steps the requested number of timesteps. */
2094     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2095     if (ts->steps >= ts->max_steps)
2096       ts->reason = TS_CONVERGED_ITS;
2097     else if (ts->ptime >= ts->max_time)
2098       ts->reason = TS_CONVERGED_TIME;
2099     while (!ts->reason) {
2100       ierr = TSStep(ts);CHKERRQ(ierr);
2101       ierr = TSPostStep(ts);CHKERRQ(ierr);
2102       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2103     }
2104     if (ts->exact_final_time && ts->ptime > ts->max_time) {
2105       ierr = TSInterpolate(ts,ts->max_time,x);CHKERRQ(ierr);
2106       if (ftime) *ftime = ts->max_time;
2107     } else {
2108       ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
2109       if (ftime) *ftime = ts->ptime;
2110     }
2111   }
2112   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2113   if (flg && !PetscPreLoadingOn) {
2114     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr);
2115     ierr = TSView(ts,viewer);CHKERRQ(ierr);
2116     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2117   }
2118   PetscFunctionReturn(0);
2119 }
2120 
2121 #undef __FUNCT__
2122 #define __FUNCT__ "TSMonitor"
2123 /*@
2124    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
2125 
2126    Collective on TS
2127 
2128    Input Parameters:
2129 +  ts - time stepping context obtained from TSCreate()
2130 .  step - step number that has just completed
2131 .  ptime - model time of the state
2132 -  x - state at the current model time
2133 
2134    Notes:
2135    TSMonitor() is typically used within the time stepping implementations.
2136    Users might call this function when using the TSStep() interface instead of TSSolve().
2137 
2138    Level: advanced
2139 
2140 .keywords: TS, timestep
2141 @*/
2142 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
2143 {
2144   PetscErrorCode ierr;
2145   PetscInt       i,n = ts->numbermonitors;
2146 
2147   PetscFunctionBegin;
2148   for (i=0; i<n; i++) {
2149     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
2150   }
2151   PetscFunctionReturn(0);
2152 }
2153 
2154 /* ------------------------------------------------------------------------*/
2155 
2156 #undef __FUNCT__
2157 #define __FUNCT__ "TSMonitorLGCreate"
2158 /*@C
2159    TSMonitorLGCreate - Creates a line graph context for use with
2160    TS to monitor convergence of preconditioned residual norms.
2161 
2162    Collective on TS
2163 
2164    Input Parameters:
2165 +  host - the X display to open, or null for the local machine
2166 .  label - the title to put in the title bar
2167 .  x, y - the screen coordinates of the upper left coordinate of the window
2168 -  m, n - the screen width and height in pixels
2169 
2170    Output Parameter:
2171 .  draw - the drawing context
2172 
2173    Options Database Key:
2174 .  -ts_monitor_draw - automatically sets line graph monitor
2175 
2176    Notes:
2177    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
2178 
2179    Level: intermediate
2180 
2181 .keywords: TS, monitor, line graph, residual, seealso
2182 
2183 .seealso: TSMonitorLGDestroy(), TSMonitorSet()
2184 
2185 @*/
2186 PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2187 {
2188   PetscDraw      win;
2189   PetscErrorCode ierr;
2190 
2191   PetscFunctionBegin;
2192   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
2193   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
2194   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
2195   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
2196 
2197   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
2198   PetscFunctionReturn(0);
2199 }
2200 
2201 #undef __FUNCT__
2202 #define __FUNCT__ "TSMonitorLG"
2203 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
2204 {
2205   PetscDrawLG    lg = (PetscDrawLG) monctx;
2206   PetscReal      x,y = ptime;
2207   PetscErrorCode ierr;
2208 
2209   PetscFunctionBegin;
2210   if (!monctx) {
2211     MPI_Comm    comm;
2212     PetscViewer viewer;
2213 
2214     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
2215     viewer = PETSC_VIEWER_DRAW_(comm);
2216     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
2217   }
2218 
2219   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2220   x = (PetscReal)n;
2221   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2222   if (n < 20 || (n % 5)) {
2223     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2224   }
2225   PetscFunctionReturn(0);
2226 }
2227 
2228 #undef __FUNCT__
2229 #define __FUNCT__ "TSMonitorLGDestroy"
2230 /*@C
2231    TSMonitorLGDestroy - Destroys a line graph context that was created
2232    with TSMonitorLGCreate().
2233 
2234    Collective on PetscDrawLG
2235 
2236    Input Parameter:
2237 .  draw - the drawing context
2238 
2239    Level: intermediate
2240 
2241 .keywords: TS, monitor, line graph, destroy
2242 
2243 .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
2244 @*/
2245 PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
2246 {
2247   PetscDraw      draw;
2248   PetscErrorCode ierr;
2249 
2250   PetscFunctionBegin;
2251   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
2252   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
2253   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
2254   PetscFunctionReturn(0);
2255 }
2256 
2257 #undef __FUNCT__
2258 #define __FUNCT__ "TSGetTime"
2259 /*@
2260    TSGetTime - Gets the time of the most recently completed step.
2261 
2262    Not Collective
2263 
2264    Input Parameter:
2265 .  ts - the TS context obtained from TSCreate()
2266 
2267    Output Parameter:
2268 .  t  - the current time
2269 
2270    Level: beginner
2271 
2272    Note:
2273    When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
2274    TSSetPreStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.
2275 
2276 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
2277 
2278 .keywords: TS, get, time
2279 @*/
2280 PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
2281 {
2282   PetscFunctionBegin;
2283   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2284   PetscValidRealPointer(t,2);
2285   *t = ts->ptime;
2286   PetscFunctionReturn(0);
2287 }
2288 
2289 #undef __FUNCT__
2290 #define __FUNCT__ "TSSetTime"
2291 /*@
2292    TSSetTime - Allows one to reset the time.
2293 
2294    Logically Collective on TS
2295 
2296    Input Parameters:
2297 +  ts - the TS context obtained from TSCreate()
2298 -  time - the time
2299 
2300    Level: intermediate
2301 
2302 .seealso: TSGetTime(), TSSetDuration()
2303 
2304 .keywords: TS, set, time
2305 @*/
2306 PetscErrorCode  TSSetTime(TS ts, PetscReal t)
2307 {
2308   PetscFunctionBegin;
2309   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2310   PetscValidLogicalCollectiveReal(ts,t,2);
2311   ts->ptime = t;
2312   PetscFunctionReturn(0);
2313 }
2314 
2315 #undef __FUNCT__
2316 #define __FUNCT__ "TSSetOptionsPrefix"
2317 /*@C
2318    TSSetOptionsPrefix - Sets the prefix used for searching for all
2319    TS options in the database.
2320 
2321    Logically Collective on TS
2322 
2323    Input Parameter:
2324 +  ts     - The TS context
2325 -  prefix - The prefix to prepend to all option names
2326 
2327    Notes:
2328    A hyphen (-) must NOT be given at the beginning of the prefix name.
2329    The first character of all runtime options is AUTOMATICALLY the
2330    hyphen.
2331 
2332    Level: advanced
2333 
2334 .keywords: TS, set, options, prefix, database
2335 
2336 .seealso: TSSetFromOptions()
2337 
2338 @*/
2339 PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2340 {
2341   PetscErrorCode ierr;
2342   SNES           snes;
2343 
2344   PetscFunctionBegin;
2345   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2346   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2347   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2348   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2349   PetscFunctionReturn(0);
2350 }
2351 
2352 
2353 #undef __FUNCT__
2354 #define __FUNCT__ "TSAppendOptionsPrefix"
2355 /*@C
2356    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2357    TS options in the database.
2358 
2359    Logically Collective on TS
2360 
2361    Input Parameter:
2362 +  ts     - The TS context
2363 -  prefix - The prefix to prepend to all option names
2364 
2365    Notes:
2366    A hyphen (-) must NOT be given at the beginning of the prefix name.
2367    The first character of all runtime options is AUTOMATICALLY the
2368    hyphen.
2369 
2370    Level: advanced
2371 
2372 .keywords: TS, append, options, prefix, database
2373 
2374 .seealso: TSGetOptionsPrefix()
2375 
2376 @*/
2377 PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2378 {
2379   PetscErrorCode ierr;
2380   SNES           snes;
2381 
2382   PetscFunctionBegin;
2383   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2384   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2385   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2386   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2387   PetscFunctionReturn(0);
2388 }
2389 
2390 #undef __FUNCT__
2391 #define __FUNCT__ "TSGetOptionsPrefix"
2392 /*@C
2393    TSGetOptionsPrefix - Sets the prefix used for searching for all
2394    TS options in the database.
2395 
2396    Not Collective
2397 
2398    Input Parameter:
2399 .  ts - The TS context
2400 
2401    Output Parameter:
2402 .  prefix - A pointer to the prefix string used
2403 
2404    Notes: On the fortran side, the user should pass in a string 'prifix' of
2405    sufficient length to hold the prefix.
2406 
2407    Level: intermediate
2408 
2409 .keywords: TS, get, options, prefix, database
2410 
2411 .seealso: TSAppendOptionsPrefix()
2412 @*/
2413 PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2414 {
2415   PetscErrorCode ierr;
2416 
2417   PetscFunctionBegin;
2418   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2419   PetscValidPointer(prefix,2);
2420   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2421   PetscFunctionReturn(0);
2422 }
2423 
2424 #undef __FUNCT__
2425 #define __FUNCT__ "TSGetRHSJacobian"
2426 /*@C
2427    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2428 
2429    Not Collective, but parallel objects are returned if TS is parallel
2430 
2431    Input Parameter:
2432 .  ts  - The TS context obtained from TSCreate()
2433 
2434    Output Parameters:
2435 +  J   - The Jacobian J of F, where U_t = F(U,t)
2436 .  M   - The preconditioner matrix, usually the same as J
2437 .  func - Function to compute the Jacobian of the RHS
2438 -  ctx - User-defined context for Jacobian evaluation routine
2439 
2440    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2441 
2442    Level: intermediate
2443 
2444 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2445 
2446 .keywords: TS, timestep, get, matrix, Jacobian
2447 @*/
2448 PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2449 {
2450   PetscErrorCode ierr;
2451   SNES           snes;
2452   DM             dm;
2453 
2454   PetscFunctionBegin;
2455   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2456   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2457   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2458   ierr = DMTSGetRHSJacobian(dm,func,ctx);CHKERRQ(ierr);
2459   PetscFunctionReturn(0);
2460 }
2461 
2462 #undef __FUNCT__
2463 #define __FUNCT__ "TSGetIJacobian"
2464 /*@C
2465    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
2466 
2467    Not Collective, but parallel objects are returned if TS is parallel
2468 
2469    Input Parameter:
2470 .  ts  - The TS context obtained from TSCreate()
2471 
2472    Output Parameters:
2473 +  A   - The Jacobian of F(t,U,U_t)
2474 .  B   - The preconditioner matrix, often the same as A
2475 .  f   - The function to compute the matrices
2476 - ctx - User-defined context for Jacobian evaluation routine
2477 
2478    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2479 
2480    Level: advanced
2481 
2482 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2483 
2484 .keywords: TS, timestep, get, matrix, Jacobian
2485 @*/
2486 PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2487 {
2488   PetscErrorCode ierr;
2489   SNES           snes;
2490   DM             dm;
2491   PetscFunctionBegin;
2492   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2493   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2494   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2495   ierr = DMTSGetIJacobian(dm,f,ctx);CHKERRQ(ierr);
2496   PetscFunctionReturn(0);
2497 }
2498 
2499 typedef struct {
2500   PetscViewer viewer;
2501   Vec         initialsolution;
2502   PetscBool   showinitial;
2503 } TSMonitorSolutionCtx;
2504 
2505 #undef __FUNCT__
2506 #define __FUNCT__ "TSMonitorSolution"
2507 /*@C
2508    TSMonitorSolution - Monitors progress of the TS solvers by calling
2509    VecView() for the solution at each timestep
2510 
2511    Collective on TS
2512 
2513    Input Parameters:
2514 +  ts - the TS context
2515 .  step - current time-step
2516 .  ptime - current time
2517 -  dummy - either a viewer or PETSC_NULL
2518 
2519    Level: intermediate
2520 
2521 .keywords: TS,  vector, monitor, view
2522 
2523 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2524 @*/
2525 PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2526 {
2527   PetscErrorCode       ierr;
2528   TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy;
2529 
2530   PetscFunctionBegin;
2531   if (!step && ictx->showinitial) {
2532     if (!ictx->initialsolution) {
2533       ierr = VecDuplicate(x,&ictx->initialsolution);CHKERRQ(ierr);
2534     }
2535     ierr = VecCopy(x,ictx->initialsolution);CHKERRQ(ierr);
2536   }
2537   if (ictx->showinitial) {
2538     PetscReal pause;
2539     ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr);
2540     ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr);
2541     ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr);
2542     ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr);
2543     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr);
2544   }
2545   ierr = VecView(x,ictx->viewer);CHKERRQ(ierr);
2546   if (ictx->showinitial) {
2547     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr);
2548   }
2549   PetscFunctionReturn(0);
2550 }
2551 
2552 
2553 #undef __FUNCT__
2554 #define __FUNCT__ "TSMonitorSolutionDestroy"
2555 /*@C
2556    TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution
2557 
2558    Collective on TS
2559 
2560    Input Parameters:
2561 .    ctx - the monitor context
2562 
2563    Level: intermediate
2564 
2565 .keywords: TS,  vector, monitor, view
2566 
2567 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
2568 @*/
2569 PetscErrorCode  TSMonitorSolutionDestroy(void **ctx)
2570 {
2571   PetscErrorCode       ierr;
2572   TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx;
2573 
2574   PetscFunctionBegin;
2575   ierr = PetscViewerDestroy(&ictx->viewer);CHKERRQ(ierr);
2576   ierr = VecDestroy(&ictx->initialsolution);CHKERRQ(ierr);
2577   ierr = PetscFree(ictx);CHKERRQ(ierr);
2578   PetscFunctionReturn(0);
2579 }
2580 
2581 #undef __FUNCT__
2582 #define __FUNCT__ "TSMonitorSolutionCreate"
2583 /*@C
2584    TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution
2585 
2586    Collective on TS
2587 
2588    Input Parameter:
2589 .    ts - time-step context
2590 
2591    Output Patameter:
2592 .    ctx - the monitor context
2593 
2594    Level: intermediate
2595 
2596 .keywords: TS,  vector, monitor, view
2597 
2598 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
2599 @*/
2600 PetscErrorCode  TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx)
2601 {
2602   PetscErrorCode       ierr;
2603   TSMonitorSolutionCtx *ictx;
2604 
2605   PetscFunctionBegin;
2606   ierr = PetscNew(TSMonitorSolutionCtx,&ictx);CHKERRQ(ierr);
2607   *ctx = (void*)ictx;
2608   if (!viewer) {
2609     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2610   }
2611   ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
2612   ictx->viewer      = viewer;
2613   ictx->showinitial = PETSC_FALSE;
2614   ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);CHKERRQ(ierr);
2615   PetscFunctionReturn(0);
2616 }
2617 
2618 #undef __FUNCT__
2619 #define __FUNCT__ "TSSetDM"
2620 /*@
2621    TSSetDM - Sets the DM that may be used by some preconditioners
2622 
2623    Logically Collective on TS and DM
2624 
2625    Input Parameters:
2626 +  ts - the preconditioner context
2627 -  dm - the dm
2628 
2629    Level: intermediate
2630 
2631 
2632 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2633 @*/
2634 PetscErrorCode  TSSetDM(TS ts,DM dm)
2635 {
2636   PetscErrorCode ierr;
2637   SNES           snes;
2638   TSDM           tsdm;
2639 
2640   PetscFunctionBegin;
2641   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2642   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
2643   if (ts->dm) {               /* Move the TSDM context over to the new DM unless the new DM already has one */
2644     PetscContainer oldcontainer,container;
2645     ierr = PetscObjectQuery((PetscObject)ts->dm,"TSDM",(PetscObject*)&oldcontainer);CHKERRQ(ierr);
2646     ierr = PetscObjectQuery((PetscObject)dm,"TSDM",(PetscObject*)&container);CHKERRQ(ierr);
2647     if (oldcontainer && !container) {
2648       ierr = DMTSCopyContext(ts->dm,dm);CHKERRQ(ierr);
2649       ierr = DMTSGetContext(ts->dm,&tsdm);CHKERRQ(ierr);
2650       if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
2651         tsdm->originaldm = dm;
2652       }
2653     }
2654     ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
2655   }
2656   ts->dm = dm;
2657   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2658   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
2659   PetscFunctionReturn(0);
2660 }
2661 
2662 #undef __FUNCT__
2663 #define __FUNCT__ "TSGetDM"
2664 /*@
2665    TSGetDM - Gets the DM that may be used by some preconditioners
2666 
2667    Not Collective
2668 
2669    Input Parameter:
2670 . ts - the preconditioner context
2671 
2672    Output Parameter:
2673 .  dm - the dm
2674 
2675    Level: intermediate
2676 
2677 
2678 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2679 @*/
2680 PetscErrorCode  TSGetDM(TS ts,DM *dm)
2681 {
2682   PetscErrorCode ierr;
2683 
2684   PetscFunctionBegin;
2685   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2686   if (!ts->dm) {
2687     ierr = DMShellCreate(((PetscObject)ts)->comm,&ts->dm);CHKERRQ(ierr);
2688     if (ts->snes) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
2689   }
2690   *dm = ts->dm;
2691   PetscFunctionReturn(0);
2692 }
2693 
2694 #undef __FUNCT__
2695 #define __FUNCT__ "SNESTSFormFunction"
2696 /*@
2697    SNESTSFormFunction - Function to evaluate nonlinear residual
2698 
2699    Logically Collective on SNES
2700 
2701    Input Parameter:
2702 + snes - nonlinear solver
2703 . X - the current state at which to evaluate the residual
2704 - ctx - user context, must be a TS
2705 
2706    Output Parameter:
2707 . F - the nonlinear residual
2708 
2709    Notes:
2710    This function is not normally called by users and is automatically registered with the SNES used by TS.
2711    It is most frequently passed to MatFDColoringSetFunction().
2712 
2713    Level: advanced
2714 
2715 .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2716 @*/
2717 PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2718 {
2719   TS ts = (TS)ctx;
2720   PetscErrorCode ierr;
2721 
2722   PetscFunctionBegin;
2723   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2724   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2725   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
2726   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
2727   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
2728   PetscFunctionReturn(0);
2729 }
2730 
2731 #undef __FUNCT__
2732 #define __FUNCT__ "SNESTSFormJacobian"
2733 /*@
2734    SNESTSFormJacobian - Function to evaluate the Jacobian
2735 
2736    Collective on SNES
2737 
2738    Input Parameter:
2739 + snes - nonlinear solver
2740 . X - the current state at which to evaluate the residual
2741 - ctx - user context, must be a TS
2742 
2743    Output Parameter:
2744 + A - the Jacobian
2745 . B - the preconditioning matrix (may be the same as A)
2746 - flag - indicates any structure change in the matrix
2747 
2748    Notes:
2749    This function is not normally called by users and is automatically registered with the SNES used by TS.
2750 
2751    Level: developer
2752 
2753 .seealso: SNESSetJacobian()
2754 @*/
2755 PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2756 {
2757   TS ts = (TS)ctx;
2758   PetscErrorCode ierr;
2759 
2760   PetscFunctionBegin;
2761   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2762   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2763   PetscValidPointer(A,3);
2764   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
2765   PetscValidPointer(B,4);
2766   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
2767   PetscValidPointer(flag,5);
2768   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
2769   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
2770   PetscFunctionReturn(0);
2771 }
2772 
2773 #undef __FUNCT__
2774 #define __FUNCT__ "TSComputeRHSFunctionLinear"
2775 /*@C
2776    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only
2777 
2778    Collective on TS
2779 
2780    Input Arguments:
2781 +  ts - time stepping context
2782 .  t - time at which to evaluate
2783 .  X - state at which to evaluate
2784 -  ctx - context
2785 
2786    Output Arguments:
2787 .  F - right hand side
2788 
2789    Level: intermediate
2790 
2791    Notes:
2792    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
2793    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
2794 
2795 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
2796 @*/
2797 PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
2798 {
2799   PetscErrorCode ierr;
2800   Mat Arhs,Brhs;
2801   MatStructure flg2;
2802 
2803   PetscFunctionBegin;
2804   ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
2805   ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
2806   ierr = MatMult(Arhs,X,F);CHKERRQ(ierr);
2807   PetscFunctionReturn(0);
2808 }
2809 
2810 #undef __FUNCT__
2811 #define __FUNCT__ "TSComputeRHSJacobianConstant"
2812 /*@C
2813    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
2814 
2815    Collective on TS
2816 
2817    Input Arguments:
2818 +  ts - time stepping context
2819 .  t - time at which to evaluate
2820 .  X - state at which to evaluate
2821 -  ctx - context
2822 
2823    Output Arguments:
2824 +  A - pointer to operator
2825 .  B - pointer to preconditioning matrix
2826 -  flg - matrix structure flag
2827 
2828    Level: intermediate
2829 
2830    Notes:
2831    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
2832 
2833 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
2834 @*/
2835 PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2836 {
2837 
2838   PetscFunctionBegin;
2839   *flg = SAME_PRECONDITIONER;
2840   PetscFunctionReturn(0);
2841 }
2842 
2843 #undef __FUNCT__
2844 #define __FUNCT__ "TSComputeIFunctionLinear"
2845 /*@C
2846    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
2847 
2848    Collective on TS
2849 
2850    Input Arguments:
2851 +  ts - time stepping context
2852 .  t - time at which to evaluate
2853 .  X - state at which to evaluate
2854 .  Xdot - time derivative of state vector
2855 -  ctx - context
2856 
2857    Output Arguments:
2858 .  F - left hand side
2859 
2860    Level: intermediate
2861 
2862    Notes:
2863    The assumption here is that the left hand side is of the form A*Xdot (and not A*Xdot + B*X). For other cases, the
2864    user is required to write their own TSComputeIFunction.
2865    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
2866    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
2867 
2868 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
2869 @*/
2870 PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
2871 {
2872   PetscErrorCode ierr;
2873   Mat A,B;
2874   MatStructure flg2;
2875 
2876   PetscFunctionBegin;
2877   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2878   ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr);
2879   ierr = MatMult(A,Xdot,F);CHKERRQ(ierr);
2880   PetscFunctionReturn(0);
2881 }
2882 
2883 #undef __FUNCT__
2884 #define __FUNCT__ "TSComputeIJacobianConstant"
2885 /*@C
2886    TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent.
2887 
2888    Collective on TS
2889 
2890    Input Arguments:
2891 +  ts - time stepping context
2892 .  t - time at which to evaluate
2893 .  X - state at which to evaluate
2894 .  Xdot - time derivative of state vector
2895 .  shift - shift to apply
2896 -  ctx - context
2897 
2898    Output Arguments:
2899 +  A - pointer to operator
2900 .  B - pointer to preconditioning matrix
2901 -  flg - matrix structure flag
2902 
2903    Level: intermediate
2904 
2905    Notes:
2906    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
2907 
2908 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
2909 @*/
2910 PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2911 {
2912 
2913   PetscFunctionBegin;
2914   *flg = SAME_PRECONDITIONER;
2915   PetscFunctionReturn(0);
2916 }
2917 
2918 
2919 #undef __FUNCT__
2920 #define __FUNCT__ "TSGetConvergedReason"
2921 /*@
2922    TSGetConvergedReason - Gets the reason the TS iteration was stopped.
2923 
2924    Not Collective
2925 
2926    Input Parameter:
2927 .  ts - the TS context
2928 
2929    Output Parameter:
2930 .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
2931             manual pages for the individual convergence tests for complete lists
2932 
2933    Level: intermediate
2934 
2935    Notes:
2936    Can only be called after the call to TSSolve() is complete.
2937 
2938 .keywords: TS, nonlinear, set, convergence, test
2939 
2940 .seealso: TSSetConvergenceTest(), TSConvergedReason
2941 @*/
2942 PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
2943 {
2944   PetscFunctionBegin;
2945   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2946   PetscValidPointer(reason,2);
2947   *reason = ts->reason;
2948   PetscFunctionReturn(0);
2949 }
2950 
2951 #undef __FUNCT__
2952 #define __FUNCT__ "TSGetSNESIterations"
2953 /*@
2954    TSGetSNESIterations - Gets the total number of nonlinear iterations
2955    used by the time integrator.
2956 
2957    Not Collective
2958 
2959    Input Parameter:
2960 .  ts - TS context
2961 
2962    Output Parameter:
2963 .  nits - number of nonlinear iterations
2964 
2965    Notes:
2966    This counter is reset to zero for each successive call to TSSolve().
2967 
2968    Level: intermediate
2969 
2970 .keywords: TS, get, number, nonlinear, iterations
2971 
2972 .seealso:  TSGetKSPIterations()
2973 @*/
2974 PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
2975 {
2976   PetscFunctionBegin;
2977   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2978   PetscValidIntPointer(nits,2);
2979   *nits = ts->snes_its;
2980   PetscFunctionReturn(0);
2981 }
2982 
2983 #undef __FUNCT__
2984 #define __FUNCT__ "TSGetKSPIterations"
2985 /*@
2986    TSGetKSPIterations - Gets the total number of linear iterations
2987    used by the time integrator.
2988 
2989    Not Collective
2990 
2991    Input Parameter:
2992 .  ts - TS context
2993 
2994    Output Parameter:
2995 .  lits - number of linear iterations
2996 
2997    Notes:
2998    This counter is reset to zero for each successive call to TSSolve().
2999 
3000    Level: intermediate
3001 
3002 .keywords: TS, get, number, linear, iterations
3003 
3004 .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
3005 @*/
3006 PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
3007 {
3008   PetscFunctionBegin;
3009   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3010   PetscValidIntPointer(lits,2);
3011   *lits = ts->ksp_its;
3012   PetscFunctionReturn(0);
3013 }
3014 
3015 #undef __FUNCT__
3016 #define __FUNCT__ "TSGetStepRejections"
3017 /*@
3018    TSGetStepRejections - Gets the total number of rejected steps.
3019 
3020    Not Collective
3021 
3022    Input Parameter:
3023 .  ts - TS context
3024 
3025    Output Parameter:
3026 .  rejects - number of steps rejected
3027 
3028    Notes:
3029    This counter is reset to zero for each successive call to TSSolve().
3030 
3031    Level: intermediate
3032 
3033 .keywords: TS, get, number
3034 
3035 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
3036 @*/
3037 PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
3038 {
3039   PetscFunctionBegin;
3040   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3041   PetscValidIntPointer(rejects,2);
3042   *rejects = ts->reject;
3043   PetscFunctionReturn(0);
3044 }
3045 
3046 #undef __FUNCT__
3047 #define __FUNCT__ "TSGetSNESFailures"
3048 /*@
3049    TSGetSNESFailures - Gets the total number of failed SNES solves
3050 
3051    Not Collective
3052 
3053    Input Parameter:
3054 .  ts - TS context
3055 
3056    Output Parameter:
3057 .  fails - number of failed nonlinear solves
3058 
3059    Notes:
3060    This counter is reset to zero for each successive call to TSSolve().
3061 
3062    Level: intermediate
3063 
3064 .keywords: TS, get, number
3065 
3066 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
3067 @*/
3068 PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
3069 {
3070   PetscFunctionBegin;
3071   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3072   PetscValidIntPointer(fails,2);
3073   *fails = ts->num_snes_failures;
3074   PetscFunctionReturn(0);
3075 }
3076 
3077 #undef __FUNCT__
3078 #define __FUNCT__ "TSSetMaxStepRejections"
3079 /*@
3080    TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails
3081 
3082    Not Collective
3083 
3084    Input Parameter:
3085 +  ts - TS context
3086 -  rejects - maximum number of rejected steps, pass -1 for unlimited
3087 
3088    Notes:
3089    The counter is reset to zero for each step
3090 
3091    Options Database Key:
3092  .  -ts_max_reject - Maximum number of step rejections before a step fails
3093 
3094    Level: intermediate
3095 
3096 .keywords: TS, set, maximum, number
3097 
3098 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3099 @*/
3100 PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
3101 {
3102   PetscFunctionBegin;
3103   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3104   ts->max_reject = rejects;
3105   PetscFunctionReturn(0);
3106 }
3107 
3108 #undef __FUNCT__
3109 #define __FUNCT__ "TSSetMaxSNESFailures"
3110 /*@
3111    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves
3112 
3113    Not Collective
3114 
3115    Input Parameter:
3116 +  ts - TS context
3117 -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited
3118 
3119    Notes:
3120    The counter is reset to zero for each successive call to TSSolve().
3121 
3122    Options Database Key:
3123  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures
3124 
3125    Level: intermediate
3126 
3127 .keywords: TS, set, maximum, number
3128 
3129 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
3130 @*/
3131 PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
3132 {
3133   PetscFunctionBegin;
3134   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3135   ts->max_snes_failures = fails;
3136   PetscFunctionReturn(0);
3137 }
3138 
3139 #undef __FUNCT__
3140 #define __FUNCT__ "TSSetErrorIfStepFails()"
3141 /*@
3142    TSSetErrorIfStepFails - Error if no step succeeds
3143 
3144    Not Collective
3145 
3146    Input Parameter:
3147 +  ts - TS context
3148 -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure
3149 
3150    Options Database Key:
3151  .  -ts_error_if_step_fails - Error if no step succeeds
3152 
3153    Level: intermediate
3154 
3155 .keywords: TS, set, error
3156 
3157 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3158 @*/
3159 PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
3160 {
3161   PetscFunctionBegin;
3162   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3163   ts->errorifstepfailed = err;
3164   PetscFunctionReturn(0);
3165 }
3166 
3167 #undef __FUNCT__
3168 #define __FUNCT__ "TSMonitorSolutionBinary"
3169 /*@C
3170    TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file
3171 
3172    Collective on TS
3173 
3174    Input Parameters:
3175 +  ts - the TS context
3176 .  step - current time-step
3177 .  ptime - current time
3178 .  x - current state
3179 -  viewer - binary viewer
3180 
3181    Level: intermediate
3182 
3183 .keywords: TS,  vector, monitor, view
3184 
3185 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3186 @*/
3187 PetscErrorCode  TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer)
3188 {
3189   PetscErrorCode       ierr;
3190   PetscViewer          v = (PetscViewer)viewer;
3191 
3192   PetscFunctionBegin;
3193   ierr = VecView(x,v);CHKERRQ(ierr);
3194   PetscFunctionReturn(0);
3195 }
3196 
3197 #undef __FUNCT__
3198 #define __FUNCT__ "TSMonitorSolutionVTK"
3199 /*@C
3200    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.
3201 
3202    Collective on TS
3203 
3204    Input Parameters:
3205 +  ts - the TS context
3206 .  step - current time-step
3207 .  ptime - current time
3208 .  x - current state
3209 -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3210 
3211    Level: intermediate
3212 
3213    Notes:
3214    The VTK format does not allow writing multiple time steps in the same file, therefore a different file will be written for each time step.
3215    These are named according to the file name template.
3216 
3217    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().
3218 
3219 .keywords: TS,  vector, monitor, view
3220 
3221 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3222 @*/
3223 PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate)
3224 {
3225   PetscErrorCode ierr;
3226   char           filename[PETSC_MAX_PATH_LEN];
3227   PetscViewer    viewer;
3228 
3229   PetscFunctionBegin;
3230   ierr = PetscSNPrintf(filename,sizeof filename,(const char*)filenametemplate,step);CHKERRQ(ierr);
3231   ierr = PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
3232   ierr = VecView(x,viewer);CHKERRQ(ierr);
3233   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3234   PetscFunctionReturn(0);
3235 }
3236 
3237 #undef __FUNCT__
3238 #define __FUNCT__ "TSMonitorSolutionVTKDestroy"
3239 /*@C
3240    TSMonitorSolutionVTKDestroy - Destroy context for monitoring
3241 
3242    Collective on TS
3243 
3244    Input Parameters:
3245 .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3246 
3247    Level: intermediate
3248 
3249    Note:
3250    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().
3251 
3252 .keywords: TS,  vector, monitor, view
3253 
3254 .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
3255 @*/
3256 PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
3257 {
3258   PetscErrorCode ierr;
3259 
3260   PetscFunctionBegin;
3261   ierr = PetscFree(*(char**)filenametemplate);CHKERRQ(ierr);
3262   PetscFunctionReturn(0);
3263 }
3264 
3265 #undef __FUNCT__
3266 #define __FUNCT__ "TSGetAdapt"
3267 /*@
3268    TSGetAdapt - Get the adaptive controller context for the current method
3269 
3270    Collective on TS if controller has not been created yet
3271 
3272    Input Arguments:
3273 .  ts - time stepping context
3274 
3275    Output Arguments:
3276 .  adapt - adaptive controller
3277 
3278    Level: intermediate
3279 
3280 .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
3281 @*/
3282 PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
3283 {
3284   PetscErrorCode ierr;
3285 
3286   PetscFunctionBegin;
3287   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3288   PetscValidPointer(adapt,2);
3289   if (!ts->adapt) {
3290     ierr = TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);CHKERRQ(ierr);
3291     ierr = PetscLogObjectParent(ts,ts->adapt);CHKERRQ(ierr);
3292     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);CHKERRQ(ierr);
3293   }
3294   *adapt = ts->adapt;
3295   PetscFunctionReturn(0);
3296 }
3297 
3298 #undef __FUNCT__
3299 #define __FUNCT__ "TSSetTolerances"
3300 /*@
3301    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller
3302 
3303    Logically Collective
3304 
3305    Input Arguments:
3306 +  ts - time integration context
3307 .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
3308 .  vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present
3309 .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
3310 -  vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present
3311 
3312    Level: beginner
3313 
3314 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
3315 @*/
3316 PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
3317 {
3318   PetscErrorCode ierr;
3319 
3320   PetscFunctionBegin;
3321   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
3322   if (vatol) {
3323     ierr = PetscObjectReference((PetscObject)vatol);CHKERRQ(ierr);
3324     ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
3325     ts->vatol = vatol;
3326   }
3327   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
3328   if (vrtol) {
3329     ierr = PetscObjectReference((PetscObject)vrtol);CHKERRQ(ierr);
3330     ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
3331     ts->vrtol = vrtol;
3332   }
3333   PetscFunctionReturn(0);
3334 }
3335 
3336 #undef __FUNCT__
3337 #define __FUNCT__ "TSGetTolerances"
3338 /*@
3339    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
3340 
3341    Logically Collective
3342 
3343    Input Arguments:
3344 .  ts - time integration context
3345 
3346    Output Arguments:
3347 +  atol - scalar absolute tolerances, PETSC_NULL to ignore
3348 .  vatol - vector of absolute tolerances, PETSC_NULL to ignore
3349 .  rtol - scalar relative tolerances, PETSC_NULL to ignore
3350 -  vrtol - vector of relative tolerances, PETSC_NULL to ignore
3351 
3352    Level: beginner
3353 
3354 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
3355 @*/
3356 PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
3357 {
3358 
3359   PetscFunctionBegin;
3360   if (atol)  *atol  = ts->atol;
3361   if (vatol) *vatol = ts->vatol;
3362   if (rtol)  *rtol  = ts->rtol;
3363   if (vrtol) *vrtol = ts->vrtol;
3364   PetscFunctionReturn(0);
3365 }
3366 
3367 #undef __FUNCT__
3368 #define __FUNCT__ "TSErrorNormWRMS"
3369 /*@
3370    TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state
3371 
3372    Collective on TS
3373 
3374    Input Arguments:
3375 +  ts - time stepping context
3376 -  Y - state vector to be compared to ts->vec_sol
3377 
3378    Output Arguments:
3379 .  norm - weighted norm, a value of 1.0 is considered small
3380 
3381    Level: developer
3382 
3383 .seealso: TSSetTolerances()
3384 @*/
3385 PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
3386 {
3387   PetscErrorCode ierr;
3388   PetscInt i,n,N;
3389   const PetscScalar *x,*y;
3390   Vec X;
3391   PetscReal sum,gsum;
3392 
3393   PetscFunctionBegin;
3394   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3395   PetscValidHeaderSpecific(Y,VEC_CLASSID,2);
3396   PetscValidPointer(norm,3);
3397   X = ts->vec_sol;
3398   PetscCheckSameTypeAndComm(X,1,Y,2);
3399   if (X == Y) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector");
3400 
3401   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
3402   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
3403   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
3404   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
3405   sum = 0.;
3406   if (ts->vatol && ts->vrtol) {
3407     const PetscScalar *atol,*rtol;
3408     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3409     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3410     for (i=0; i<n; i++) {
3411       PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3412       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3413     }
3414     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3415     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3416   } else if (ts->vatol) {       /* vector atol, scalar rtol */
3417     const PetscScalar *atol;
3418     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3419     for (i=0; i<n; i++) {
3420       PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3421       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3422     }
3423     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3424   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
3425     const PetscScalar *rtol;
3426     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3427     for (i=0; i<n; i++) {
3428       PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3429       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3430     }
3431     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3432   } else {                      /* scalar atol, scalar rtol */
3433     for (i=0; i<n; i++) {
3434       PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3435       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3436     }
3437   }
3438   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
3439   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
3440 
3441   ierr = MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);CHKERRQ(ierr);
3442   *norm = PetscSqrtReal(gsum / N);
3443   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
3444   PetscFunctionReturn(0);
3445 }
3446 
3447 #undef __FUNCT__
3448 #define __FUNCT__ "TSSetCFLTimeLocal"
3449 /*@
3450    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
3451 
3452    Logically Collective on TS
3453 
3454    Input Arguments:
3455 +  ts - time stepping context
3456 -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)
3457 
3458    Note:
3459    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
3460 
3461    Level: intermediate
3462 
3463 .seealso: TSGetCFLTime(), TSADAPTCFL
3464 @*/
3465 PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
3466 {
3467 
3468   PetscFunctionBegin;
3469   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3470   ts->cfltime_local = cfltime;
3471   ts->cfltime = -1.;
3472   PetscFunctionReturn(0);
3473 }
3474 
3475 #undef __FUNCT__
3476 #define __FUNCT__ "TSGetCFLTime"
3477 /*@
3478    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
3479 
3480    Collective on TS
3481 
3482    Input Arguments:
3483 .  ts - time stepping context
3484 
3485    Output Arguments:
3486 .  cfltime - maximum stable time step for forward Euler
3487 
3488    Level: advanced
3489 
3490 .seealso: TSSetCFLTimeLocal()
3491 @*/
3492 PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
3493 {
3494   PetscErrorCode ierr;
3495 
3496   PetscFunctionBegin;
3497   if (ts->cfltime < 0) {
3498     ierr = MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr);
3499   }
3500   *cfltime = ts->cfltime;
3501   PetscFunctionReturn(0);
3502 }
3503 
3504 #undef __FUNCT__
3505 #define __FUNCT__ "TSVISetVariableBounds"
3506 /*@
3507    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
3508 
3509    Input Parameters:
3510 .  ts   - the TS context.
3511 .  xl   - lower bound.
3512 .  xu   - upper bound.
3513 
3514    Notes:
3515    If this routine is not called then the lower and upper bounds are set to
3516    SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp().
3517 
3518    Level: advanced
3519 
3520 @*/
3521 PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
3522 {
3523   PetscErrorCode ierr;
3524   SNES           snes;
3525 
3526   PetscFunctionBegin;
3527   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3528   ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr);
3529   PetscFunctionReturn(0);
3530 }
3531 
3532 #if defined(PETSC_HAVE_MATLAB_ENGINE)
3533 #include <mex.h>
3534 
3535 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
3536 
3537 #undef __FUNCT__
3538 #define __FUNCT__ "TSComputeFunction_Matlab"
3539 /*
3540    TSComputeFunction_Matlab - Calls the function that has been set with
3541                          TSSetFunctionMatlab().
3542 
3543    Collective on TS
3544 
3545    Input Parameters:
3546 +  snes - the TS context
3547 -  x - input vector
3548 
3549    Output Parameter:
3550 .  y - function vector, as set by TSSetFunction()
3551 
3552    Notes:
3553    TSComputeFunction() is typically used within nonlinear solvers
3554    implementations, so most users would not generally call this routine
3555    themselves.
3556 
3557    Level: developer
3558 
3559 .keywords: TS, nonlinear, compute, function
3560 
3561 .seealso: TSSetFunction(), TSGetFunction()
3562 */
3563 PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
3564 {
3565   PetscErrorCode   ierr;
3566   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3567   int              nlhs = 1,nrhs = 7;
3568   mxArray          *plhs[1],*prhs[7];
3569   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
3570 
3571   PetscFunctionBegin;
3572   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
3573   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3574   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
3575   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
3576   PetscCheckSameComm(snes,1,x,3);
3577   PetscCheckSameComm(snes,1,y,5);
3578 
3579   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3580   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3581   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
3582   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
3583   prhs[0] =  mxCreateDoubleScalar((double)ls);
3584   prhs[1] =  mxCreateDoubleScalar(time);
3585   prhs[2] =  mxCreateDoubleScalar((double)lx);
3586   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3587   prhs[4] =  mxCreateDoubleScalar((double)ly);
3588   prhs[5] =  mxCreateString(sctx->funcname);
3589   prhs[6] =  sctx->ctx;
3590   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
3591   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3592   mxDestroyArray(prhs[0]);
3593   mxDestroyArray(prhs[1]);
3594   mxDestroyArray(prhs[2]);
3595   mxDestroyArray(prhs[3]);
3596   mxDestroyArray(prhs[4]);
3597   mxDestroyArray(prhs[5]);
3598   mxDestroyArray(plhs[0]);
3599   PetscFunctionReturn(0);
3600 }
3601 
3602 
3603 #undef __FUNCT__
3604 #define __FUNCT__ "TSSetFunctionMatlab"
3605 /*
3606    TSSetFunctionMatlab - Sets the function evaluation routine and function
3607    vector for use by the TS routines in solving ODEs
3608    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
3609 
3610    Logically Collective on TS
3611 
3612    Input Parameters:
3613 +  ts - the TS context
3614 -  func - function evaluation routine
3615 
3616    Calling sequence of func:
3617 $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
3618 
3619    Level: beginner
3620 
3621 .keywords: TS, nonlinear, set, function
3622 
3623 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3624 */
3625 PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
3626 {
3627   PetscErrorCode  ierr;
3628   TSMatlabContext *sctx;
3629 
3630   PetscFunctionBegin;
3631   /* currently sctx is memory bleed */
3632   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3633   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3634   /*
3635      This should work, but it doesn't
3636   sctx->ctx = ctx;
3637   mexMakeArrayPersistent(sctx->ctx);
3638   */
3639   sctx->ctx = mxDuplicateArray(ctx);
3640   ierr = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
3641   PetscFunctionReturn(0);
3642 }
3643 
3644 #undef __FUNCT__
3645 #define __FUNCT__ "TSComputeJacobian_Matlab"
3646 /*
3647    TSComputeJacobian_Matlab - Calls the function that has been set with
3648                          TSSetJacobianMatlab().
3649 
3650    Collective on TS
3651 
3652    Input Parameters:
3653 +  ts - the TS context
3654 .  x - input vector
3655 .  A, B - the matrices
3656 -  ctx - user context
3657 
3658    Output Parameter:
3659 .  flag - structure of the matrix
3660 
3661    Level: developer
3662 
3663 .keywords: TS, nonlinear, compute, function
3664 
3665 .seealso: TSSetFunction(), TSGetFunction()
3666 @*/
3667 PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
3668 {
3669   PetscErrorCode  ierr;
3670   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3671   int             nlhs = 2,nrhs = 9;
3672   mxArray         *plhs[2],*prhs[9];
3673   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
3674 
3675   PetscFunctionBegin;
3676   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3677   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3678 
3679   /* call Matlab function in ctx with arguments x and y */
3680 
3681   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
3682   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3683   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
3684   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
3685   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
3686   prhs[0] =  mxCreateDoubleScalar((double)ls);
3687   prhs[1] =  mxCreateDoubleScalar((double)time);
3688   prhs[2] =  mxCreateDoubleScalar((double)lx);
3689   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3690   prhs[4] =  mxCreateDoubleScalar((double)shift);
3691   prhs[5] =  mxCreateDoubleScalar((double)lA);
3692   prhs[6] =  mxCreateDoubleScalar((double)lB);
3693   prhs[7] =  mxCreateString(sctx->funcname);
3694   prhs[8] =  sctx->ctx;
3695   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
3696   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3697   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
3698   mxDestroyArray(prhs[0]);
3699   mxDestroyArray(prhs[1]);
3700   mxDestroyArray(prhs[2]);
3701   mxDestroyArray(prhs[3]);
3702   mxDestroyArray(prhs[4]);
3703   mxDestroyArray(prhs[5]);
3704   mxDestroyArray(prhs[6]);
3705   mxDestroyArray(prhs[7]);
3706   mxDestroyArray(plhs[0]);
3707   mxDestroyArray(plhs[1]);
3708   PetscFunctionReturn(0);
3709 }
3710 
3711 
3712 #undef __FUNCT__
3713 #define __FUNCT__ "TSSetJacobianMatlab"
3714 /*
3715    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
3716    vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function
3717 
3718    Logically Collective on TS
3719 
3720    Input Parameters:
3721 +  ts - the TS context
3722 .  A,B - Jacobian matrices
3723 .  func - function evaluation routine
3724 -  ctx - user context
3725 
3726    Calling sequence of func:
3727 $    flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
3728 
3729 
3730    Level: developer
3731 
3732 .keywords: TS, nonlinear, set, function
3733 
3734 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3735 */
3736 PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
3737 {
3738   PetscErrorCode    ierr;
3739   TSMatlabContext *sctx;
3740 
3741   PetscFunctionBegin;
3742   /* currently sctx is memory bleed */
3743   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3744   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3745   /*
3746      This should work, but it doesn't
3747   sctx->ctx = ctx;
3748   mexMakeArrayPersistent(sctx->ctx);
3749   */
3750   sctx->ctx = mxDuplicateArray(ctx);
3751   ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
3752   PetscFunctionReturn(0);
3753 }
3754 
3755 #undef __FUNCT__
3756 #define __FUNCT__ "TSMonitor_Matlab"
3757 /*
3758    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
3759 
3760    Collective on TS
3761 
3762 .seealso: TSSetFunction(), TSGetFunction()
3763 @*/
3764 PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
3765 {
3766   PetscErrorCode  ierr;
3767   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3768   int             nlhs = 1,nrhs = 6;
3769   mxArray         *plhs[1],*prhs[6];
3770   long long int   lx = 0,ls = 0;
3771 
3772   PetscFunctionBegin;
3773   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3774   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
3775 
3776   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
3777   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3778   prhs[0] =  mxCreateDoubleScalar((double)ls);
3779   prhs[1] =  mxCreateDoubleScalar((double)it);
3780   prhs[2] =  mxCreateDoubleScalar((double)time);
3781   prhs[3] =  mxCreateDoubleScalar((double)lx);
3782   prhs[4] =  mxCreateString(sctx->funcname);
3783   prhs[5] =  sctx->ctx;
3784   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
3785   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3786   mxDestroyArray(prhs[0]);
3787   mxDestroyArray(prhs[1]);
3788   mxDestroyArray(prhs[2]);
3789   mxDestroyArray(prhs[3]);
3790   mxDestroyArray(prhs[4]);
3791   mxDestroyArray(plhs[0]);
3792   PetscFunctionReturn(0);
3793 }
3794 
3795 
3796 #undef __FUNCT__
3797 #define __FUNCT__ "TSMonitorSetMatlab"
3798 /*
3799    TSMonitorSetMatlab - Sets the monitor function from Matlab
3800 
3801    Level: developer
3802 
3803 .keywords: TS, nonlinear, set, function
3804 
3805 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3806 */
3807 PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
3808 {
3809   PetscErrorCode    ierr;
3810   TSMatlabContext *sctx;
3811 
3812   PetscFunctionBegin;
3813   /* currently sctx is memory bleed */
3814   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3815   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3816   /*
3817      This should work, but it doesn't
3818   sctx->ctx = ctx;
3819   mexMakeArrayPersistent(sctx->ctx);
3820   */
3821   sctx->ctx = mxDuplicateArray(ctx);
3822   ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
3823   PetscFunctionReturn(0);
3824 }
3825 #endif
3826