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