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