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