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