xref: /petsc/src/ts/interface/ts.c (revision f33bbcb64dc56389f04a9c62c8775beea4c588a2)
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   PetscReal      dt;
74   PetscBool      opt,flg;
75   PetscErrorCode ierr;
76   PetscViewer    monviewer;
77   char           monfilename[PETSC_MAX_PATH_LEN];
78   SNES           snes;
79 
80   PetscFunctionBegin;
81   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
82   ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr);
83 
84     /* Handle generic TS options */
85     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
86     ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
87     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
88     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
89     if (opt) {
90       ts->initial_time_step = ts->time_step = dt;
91     }
92     ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr);
93     ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections","",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr);
94     ierr = PetscOptionsBool("-ts_error_if_step_failed","Error if no step succeeds","",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr);
95 
96     /* Monitor options */
97     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
98     if (flg) {
99       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr);
100       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
101     }
102     opt  = PETSC_FALSE;
103     ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
104     if (opt) {
105       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
106     }
107     opt  = PETSC_FALSE;
108     ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
109     if (opt) {
110       ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
111     }
112 
113     /* Handle TS type options */
114     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
115 
116     /* Handle specific TS options */
117     if (ts->ops->setfromoptions) {
118       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
119     }
120 
121     /* process any options handlers added with PetscObjectAddOptionsHandler() */
122     ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
123   ierr = PetscOptionsEnd();CHKERRQ(ierr);
124 
125   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
126   /* Handle subobject options */
127   if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
128   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
129   PetscFunctionReturn(0);
130 }
131 
132 #undef  __FUNCT__
133 #define __FUNCT__ "TSViewFromOptions"
134 /*@
135   TSViewFromOptions - This function visualizes the ts based upon user options.
136 
137   Collective on TS
138 
139   Input Parameter:
140 . ts - The ts
141 
142   Level: intermediate
143 
144 .keywords: TS, view, options, database
145 .seealso: TSSetFromOptions(), TSView()
146 @*/
147 PetscErrorCode  TSViewFromOptions(TS ts,const char title[])
148 {
149   PetscViewer    viewer;
150   PetscDraw      draw;
151   PetscBool      opt = PETSC_FALSE;
152   char           fileName[PETSC_MAX_PATH_LEN];
153   PetscErrorCode ierr;
154 
155   PetscFunctionBegin;
156   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr);
157   if (opt && !PetscPreLoadingOn) {
158     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr);
159     ierr = TSView(ts, viewer);CHKERRQ(ierr);
160     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
161   }
162   opt = PETSC_FALSE;
163   ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr);
164   if (opt) {
165     ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
166     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
167     if (title) {
168       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
169     } else {
170       ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr);
171       ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr);
172     }
173     ierr = TSView(ts, viewer);CHKERRQ(ierr);
174     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
175     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
176     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
177   }
178   PetscFunctionReturn(0);
179 }
180 
181 #undef __FUNCT__
182 #define __FUNCT__ "TSComputeRHSJacobian"
183 /*@
184    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
185       set with TSSetRHSJacobian().
186 
187    Collective on TS and Vec
188 
189    Input Parameters:
190 +  ts - the TS context
191 .  t - current timestep
192 -  x - input vector
193 
194    Output Parameters:
195 +  A - Jacobian matrix
196 .  B - optional preconditioning matrix
197 -  flag - flag indicating matrix structure
198 
199    Notes:
200    Most users should not need to explicitly call this routine, as it
201    is used internally within the nonlinear solvers.
202 
203    See KSPSetOperators() for important information about setting the
204    flag parameter.
205 
206    Level: developer
207 
208 .keywords: SNES, compute, Jacobian, matrix
209 
210 .seealso:  TSSetRHSJacobian(), KSPSetOperators()
211 @*/
212 PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
213 {
214   PetscErrorCode ierr;
215 
216   PetscFunctionBegin;
217   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
218   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
219   PetscCheckSameComm(ts,1,X,3);
220   if (ts->ops->rhsjacobian) {
221     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
222     *flg = DIFFERENT_NONZERO_PATTERN;
223     PetscStackPush("TS user Jacobian function");
224     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
225     PetscStackPop;
226     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
227     /* make sure user returned a correct Jacobian and preconditioner */
228     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
229     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
230   } else {
231     ierr = MatZeroEntries(*A);CHKERRQ(ierr);
232     if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);}
233     *flg = SAME_NONZERO_PATTERN;
234   }
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "TSComputeRHSFunction"
240 /*@
241    TSComputeRHSFunction - Evaluates the right-hand-side function.
242 
243    Collective on TS and Vec
244 
245    Input Parameters:
246 +  ts - the TS context
247 .  t - current time
248 -  x - state vector
249 
250    Output Parameter:
251 .  y - right hand side
252 
253    Note:
254    Most users should not need to explicitly call this routine, as it
255    is used internally within the nonlinear solvers.
256 
257    Level: developer
258 
259 .keywords: TS, compute
260 
261 .seealso: TSSetRHSFunction(), TSComputeIFunction()
262 @*/
263 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
264 {
265   PetscErrorCode ierr;
266 
267   PetscFunctionBegin;
268   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
269   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
270   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
271 
272   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
273   if (ts->ops->rhsfunction) {
274     PetscStackPush("TS user right-hand-side function");
275     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
276     PetscStackPop;
277   } else {
278     ierr = VecZeroEntries(y);CHKERRQ(ierr);
279   }
280   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
281   PetscFunctionReturn(0);
282 }
283 
284 #undef __FUNCT__
285 #define __FUNCT__ "TSGetRHSVec_Private"
286 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
287 {
288   PetscErrorCode ierr;
289 
290   PetscFunctionBegin;
291   if (!ts->Frhs) {
292     ierr = VecDuplicate(ts->vec_sol,&ts->Frhs);CHKERRQ(ierr);
293   }
294   *Frhs = ts->Frhs;
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "TSGetRHSMats_Private"
300 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
301 {
302   PetscErrorCode ierr;
303   Mat A,B;
304 
305   PetscFunctionBegin;
306   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
307   if (Arhs) {
308     if (!ts->Arhs) {
309       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
310     }
311     *Arhs = ts->Arhs;
312   }
313   if (Brhs) {
314     if (!ts->Brhs) {
315       ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
316     }
317     *Brhs = ts->Brhs;
318   }
319   PetscFunctionReturn(0);
320 }
321 
322 #undef __FUNCT__
323 #define __FUNCT__ "TSComputeIFunction"
324 /*@
325    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
326 
327    Collective on TS and Vec
328 
329    Input Parameters:
330 +  ts - the TS context
331 .  t - current time
332 .  X - state vector
333 .  Xdot - time derivative of state vector
334 -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
335 
336    Output Parameter:
337 .  Y - right hand side
338 
339    Note:
340    Most users should not need to explicitly call this routine, as it
341    is used internally within the nonlinear solvers.
342 
343    If the user did did not write their equations in implicit form, this
344    function recasts them in implicit form.
345 
346    Level: developer
347 
348 .keywords: TS, compute
349 
350 .seealso: TSSetIFunction(), TSComputeRHSFunction()
351 @*/
352 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
353 {
354   PetscErrorCode ierr;
355 
356   PetscFunctionBegin;
357   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
358   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
359   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
360   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
361 
362   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
363   if (ts->ops->ifunction) {
364     PetscStackPush("TS user implicit function");
365     ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr);
366     PetscStackPop;
367   }
368   if (imex) {
369     if (!ts->ops->ifunction) {ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);}
370   } else {
371     if (!ts->ops->ifunction) {
372       ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr);
373       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
374     } else {
375       Vec Frhs;
376       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
377       ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr);
378       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
379     }
380   }
381   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
382   PetscFunctionReturn(0);
383 }
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "TSComputeIJacobian"
387 /*@
388    TSComputeIJacobian - Evaluates the Jacobian of the DAE
389 
390    Collective on TS and Vec
391 
392    Input
393       Input Parameters:
394 +  ts - the TS context
395 .  t - current timestep
396 .  X - state vector
397 .  Xdot - time derivative of state vector
398 .  shift - shift to apply, see note below
399 -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
400 
401    Output Parameters:
402 +  A - Jacobian matrix
403 .  B - optional preconditioning matrix
404 -  flag - flag indicating matrix structure
405 
406    Notes:
407    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
408 
409    dF/dX + shift*dF/dXdot
410 
411    Most users should not need to explicitly call this routine, as it
412    is used internally within the nonlinear solvers.
413 
414    Level: developer
415 
416 .keywords: TS, compute, Jacobian, matrix
417 
418 .seealso:  TSSetIJacobian()
419 @*/
420 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
421 {
422   PetscErrorCode ierr;
423 
424   PetscFunctionBegin;
425   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
426   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
427   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
428   PetscValidPointer(A,6);
429   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
430   PetscValidPointer(B,7);
431   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
432   PetscValidPointer(flg,8);
433 
434   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
435   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
436   if (ts->ops->ijacobian) {
437     PetscStackPush("TS user implicit Jacobian");
438     ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr);
439     PetscStackPop;
440     /* make sure user returned a correct Jacobian and preconditioner */
441     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
442     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
443   }
444   if (imex) {
445     if (!ts->ops->ijacobian) {  /* system was written as Xdot = F(t,X) */
446       ierr = MatZeroEntries(*A);CHKERRQ(ierr);
447       ierr = MatShift(*A,1.0);CHKERRQ(ierr);
448       if (*A != *B) {
449         ierr = MatZeroEntries(*B);CHKERRQ(ierr);
450         ierr = MatShift(*B,shift);CHKERRQ(ierr);
451       }
452       *flg = SAME_PRECONDITIONER;
453     }
454   } else {
455     if (!ts->ops->ijacobian) {
456       ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr);
457       ierr = MatScale(*A,-1);CHKERRQ(ierr);
458       ierr = MatShift(*A,shift);CHKERRQ(ierr);
459       if (*A != *B) {
460         ierr = MatScale(*B,-1);CHKERRQ(ierr);
461         ierr = MatShift(*B,shift);CHKERRQ(ierr);
462       }
463     } else if (ts->ops->rhsjacobian) {
464       Mat Arhs,Brhs;
465       MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
466       ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
467       ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
468       axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
469       ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr);
470       if (*A != *B) {
471         ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr);
472       }
473       *flg = PetscMin(*flg,flg2);
474     }
475   }
476   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
477   PetscFunctionReturn(0);
478 }
479 
480 #undef __FUNCT__
481 #define __FUNCT__ "TSSetRHSFunction"
482 /*@C
483     TSSetRHSFunction - Sets the routine for evaluating the function,
484     F(t,u), where U_t = F(t,u).
485 
486     Logically Collective on TS
487 
488     Input Parameters:
489 +   ts - the TS context obtained from TSCreate()
490 .   r - vector to put the computed right hand side (or PETSC_NULL to have it created)
491 .   f - routine for evaluating the right-hand-side function
492 -   ctx - [optional] user-defined context for private data for the
493           function evaluation routine (may be PETSC_NULL)
494 
495     Calling sequence of func:
496 $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
497 
498 +   t - current timestep
499 .   u - input vector
500 .   F - function vector
501 -   ctx - [optional] user-defined function context
502 
503     Important:
504     The user MUST call either this routine or TSSetMatrices().
505 
506     Level: beginner
507 
508 .keywords: TS, timestep, set, right-hand-side, function
509 
510 .seealso: TSSetMatrices()
511 @*/
512 PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
513 {
514   PetscErrorCode ierr;
515   SNES           snes;
516 
517   PetscFunctionBegin;
518   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
519   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
520   if (f)   ts->ops->rhsfunction = f;
521   if (ctx) ts->funP             = ctx;
522   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
523   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
524   PetscFunctionReturn(0);
525 }
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "TSSetRHSJacobian"
529 /*@C
530    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
531    where U_t = F(U,t), as well as the location to store the matrix.
532    Use TSSetMatrices() for linear problems.
533 
534    Logically Collective on TS
535 
536    Input Parameters:
537 +  ts  - the TS context obtained from TSCreate()
538 .  A   - Jacobian matrix
539 .  B   - preconditioner matrix (usually same as A)
540 .  f   - the Jacobian evaluation routine
541 -  ctx - [optional] user-defined context for private data for the
542          Jacobian evaluation routine (may be PETSC_NULL)
543 
544    Calling sequence of func:
545 $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
546 
547 +  t - current timestep
548 .  u - input vector
549 .  A - matrix A, where U_t = A(t)u
550 .  B - preconditioner matrix, usually the same as A
551 .  flag - flag indicating information about the preconditioner matrix
552           structure (same as flag in KSPSetOperators())
553 -  ctx - [optional] user-defined context for matrix evaluation routine
554 
555    Notes:
556    See KSPSetOperators() for important information about setting the flag
557    output parameter in the routine func().  Be sure to read this information!
558 
559    The routine func() takes Mat * as the matrix arguments rather than Mat.
560    This allows the matrix evaluation routine to replace A and/or B with a
561    completely new matrix structure (not just different matrix elements)
562    when appropriate, for instance, if the nonzero structure is changing
563    throughout the global iterations.
564 
565    Level: beginner
566 
567 .keywords: TS, timestep, set, right-hand-side, Jacobian
568 
569 .seealso: TSDefaultComputeJacobianColor(),
570           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices()
571 
572 @*/
573 PetscErrorCode  TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
574 {
575   PetscErrorCode ierr;
576   SNES           snes;
577 
578   PetscFunctionBegin;
579   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
580   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
581   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
582   if (A) PetscCheckSameComm(ts,1,A,2);
583   if (B) PetscCheckSameComm(ts,1,B,3);
584 
585   if (f)   ts->ops->rhsjacobian = f;
586   if (ctx) ts->jacP             = ctx;
587   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
588   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
589   PetscFunctionReturn(0);
590 }
591 
592 
593 #undef __FUNCT__
594 #define __FUNCT__ "TSSetIFunction"
595 /*@C
596    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
597 
598    Logically Collective on TS
599 
600    Input Parameters:
601 +  ts  - the TS context obtained from TSCreate()
602 .  r   - vector to hold the residual (or PETSC_NULL to have it created internally)
603 .  f   - the function evaluation routine
604 -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
605 
606    Calling sequence of f:
607 $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
608 
609 +  t   - time at step/stage being solved
610 .  u   - state vector
611 .  u_t - time derivative of state vector
612 .  F   - function vector
613 -  ctx - [optional] user-defined context for matrix evaluation routine
614 
615    Important:
616    The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices().  This routine must be used when not solving an ODE.
617 
618    Level: beginner
619 
620 .keywords: TS, timestep, set, DAE, Jacobian
621 
622 .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian()
623 @*/
624 PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
625 {
626   PetscErrorCode ierr;
627   SNES           snes;
628 
629   PetscFunctionBegin;
630   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
631   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
632   if (f)   ts->ops->ifunction = f;
633   if (ctx) ts->funP           = ctx;
634   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
635   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
636   PetscFunctionReturn(0);
637 }
638 
639 #undef __FUNCT__
640 #define __FUNCT__ "TSGetIFunction"
641 /*@C
642    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
643 
644    Not Collective
645 
646    Input Parameter:
647 .  ts - the TS context
648 
649    Output Parameter:
650 +  r - vector to hold residual (or PETSC_NULL)
651 .  func - the function to compute residual (or PETSC_NULL)
652 -  ctx - the function context (or PETSC_NULL)
653 
654    Level: advanced
655 
656 .keywords: TS, nonlinear, get, function
657 
658 .seealso: TSSetIFunction(), SNESGetFunction()
659 @*/
660 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
661 {
662   PetscErrorCode ierr;
663   SNES snes;
664 
665   PetscFunctionBegin;
666   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
667   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
668   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
669   if (func) *func = ts->ops->ifunction;
670   if (ctx)  *ctx  = ts->funP;
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNCT__
675 #define __FUNCT__ "TSGetRHSFunction"
676 /*@C
677    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
678 
679    Not Collective
680 
681    Input Parameter:
682 .  ts - the TS context
683 
684    Output Parameter:
685 +  r - vector to hold computed right hand side (or PETSC_NULL)
686 .  func - the function to compute right hand side (or PETSC_NULL)
687 -  ctx - the function context (or PETSC_NULL)
688 
689    Level: advanced
690 
691 .keywords: TS, nonlinear, get, function
692 
693 .seealso: TSSetRhsfunction(), SNESGetFunction()
694 @*/
695 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
696 {
697   PetscErrorCode ierr;
698   SNES snes;
699 
700   PetscFunctionBegin;
701   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
702   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
703   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
704   if (func) *func = ts->ops->rhsfunction;
705   if (ctx)  *ctx  = ts->funP;
706   PetscFunctionReturn(0);
707 }
708 
709 #undef __FUNCT__
710 #define __FUNCT__ "TSSetIJacobian"
711 /*@C
712    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
713         you provided with TSSetIFunction().
714 
715    Logically Collective on TS
716 
717    Input Parameters:
718 +  ts  - the TS context obtained from TSCreate()
719 .  A   - Jacobian matrix
720 .  B   - preconditioning matrix for A (may be same as A)
721 .  f   - the Jacobian evaluation routine
722 -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
723 
724    Calling sequence of f:
725 $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
726 
727 +  t    - time at step/stage being solved
728 .  U    - state vector
729 .  U_t  - time derivative of state vector
730 .  a    - shift
731 .  A    - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
732 .  B    - preconditioning matrix for A, may be same as A
733 .  flag - flag indicating information about the preconditioner matrix
734           structure (same as flag in KSPSetOperators())
735 -  ctx  - [optional] user-defined context for matrix evaluation routine
736 
737    Notes:
738    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
739 
740    The matrix dF/dU + a*dF/dU_t you provide turns out to be
741    the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
742    The time integrator internally approximates U_t by W+a*U where the positive "shift"
743    a and vector W depend on the integration method, step size, and past states. For example with
744    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
745    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
746 
747    Level: beginner
748 
749 .keywords: TS, timestep, DAE, Jacobian
750 
751 .seealso: TSSetIFunction(), TSSetRHSJacobian()
752 
753 @*/
754 PetscErrorCode  TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
755 {
756   PetscErrorCode ierr;
757   SNES           snes;
758 
759   PetscFunctionBegin;
760   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
761   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
762   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
763   if (A) PetscCheckSameComm(ts,1,A,2);
764   if (B) PetscCheckSameComm(ts,1,B,3);
765   if (f)   ts->ops->ijacobian = f;
766   if (ctx) ts->jacP           = ctx;
767   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
768   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
769   PetscFunctionReturn(0);
770 }
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "TSView"
774 /*@C
775     TSView - Prints the TS data structure.
776 
777     Collective on TS
778 
779     Input Parameters:
780 +   ts - the TS context obtained from TSCreate()
781 -   viewer - visualization context
782 
783     Options Database Key:
784 .   -ts_view - calls TSView() at end of TSStep()
785 
786     Notes:
787     The available visualization contexts include
788 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
789 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
790          output where only the first processor opens
791          the file.  All other processors send their
792          data to the first processor to print.
793 
794     The user can open an alternative visualization context with
795     PetscViewerASCIIOpen() - output to a specified file.
796 
797     Level: beginner
798 
799 .keywords: TS, timestep, view
800 
801 .seealso: PetscViewerASCIIOpen()
802 @*/
803 PetscErrorCode  TSView(TS ts,PetscViewer viewer)
804 {
805   PetscErrorCode ierr;
806   const TSType   type;
807   PetscBool      iascii,isstring,isundials;
808 
809   PetscFunctionBegin;
810   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
811   if (!viewer) {
812     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
813   }
814   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
815   PetscCheckSameComm(ts,1,viewer,2);
816 
817   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
818   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
819   if (iascii) {
820     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr);
821     if (ts->ops->view) {
822       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
823       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
824       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
825     }
826     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
827     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
828     if (ts->problem_type == TS_NONLINEAR) {
829       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
830       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->max_snes_failures);CHKERRQ(ierr);
831     }
832     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
833     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
834   } else if (isstring) {
835     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
836     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
837   }
838   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
839   ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
840   if (!isundials && ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
841   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
842   PetscFunctionReturn(0);
843 }
844 
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "TSSetApplicationContext"
848 /*@
849    TSSetApplicationContext - Sets an optional user-defined context for
850    the timesteppers.
851 
852    Logically Collective on TS
853 
854    Input Parameters:
855 +  ts - the TS context obtained from TSCreate()
856 -  usrP - optional user context
857 
858    Level: intermediate
859 
860 .keywords: TS, timestep, set, application, context
861 
862 .seealso: TSGetApplicationContext()
863 @*/
864 PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
865 {
866   PetscFunctionBegin;
867   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
868   ts->user = usrP;
869   PetscFunctionReturn(0);
870 }
871 
872 #undef __FUNCT__
873 #define __FUNCT__ "TSGetApplicationContext"
874 /*@
875     TSGetApplicationContext - Gets the user-defined context for the
876     timestepper.
877 
878     Not Collective
879 
880     Input Parameter:
881 .   ts - the TS context obtained from TSCreate()
882 
883     Output Parameter:
884 .   usrP - user context
885 
886     Level: intermediate
887 
888 .keywords: TS, timestep, get, application, context
889 
890 .seealso: TSSetApplicationContext()
891 @*/
892 PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
893 {
894   PetscFunctionBegin;
895   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
896   *(void**)usrP = ts->user;
897   PetscFunctionReturn(0);
898 }
899 
900 #undef __FUNCT__
901 #define __FUNCT__ "TSGetTimeStepNumber"
902 /*@
903    TSGetTimeStepNumber - Gets the current number of timesteps.
904 
905    Not Collective
906 
907    Input Parameter:
908 .  ts - the TS context obtained from TSCreate()
909 
910    Output Parameter:
911 .  iter - number steps so far
912 
913    Level: intermediate
914 
915 .keywords: TS, timestep, get, iteration, number
916 @*/
917 PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
918 {
919   PetscFunctionBegin;
920   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
921   PetscValidIntPointer(iter,2);
922   *iter = ts->steps;
923   PetscFunctionReturn(0);
924 }
925 
926 #undef __FUNCT__
927 #define __FUNCT__ "TSSetInitialTimeStep"
928 /*@
929    TSSetInitialTimeStep - Sets the initial timestep to be used,
930    as well as the initial time.
931 
932    Logically Collective on TS
933 
934    Input Parameters:
935 +  ts - the TS context obtained from TSCreate()
936 .  initial_time - the initial time
937 -  time_step - the size of the timestep
938 
939    Level: intermediate
940 
941 .seealso: TSSetTimeStep(), TSGetTimeStep()
942 
943 .keywords: TS, set, initial, timestep
944 @*/
945 PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
946 {
947   PetscErrorCode ierr;
948 
949   PetscFunctionBegin;
950   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
951   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
952   ts->initial_time_step = time_step;
953   ts->ptime             = initial_time;
954   PetscFunctionReturn(0);
955 }
956 
957 #undef __FUNCT__
958 #define __FUNCT__ "TSSetTimeStep"
959 /*@
960    TSSetTimeStep - Allows one to reset the timestep at any time,
961    useful for simple pseudo-timestepping codes.
962 
963    Logically Collective on TS
964 
965    Input Parameters:
966 +  ts - the TS context obtained from TSCreate()
967 -  time_step - the size of the timestep
968 
969    Level: intermediate
970 
971 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
972 
973 .keywords: TS, set, timestep
974 @*/
975 PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
976 {
977   PetscFunctionBegin;
978   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
979   PetscValidLogicalCollectiveReal(ts,time_step,2);
980   ts->time_step      = time_step;
981   ts->next_time_step = time_step;
982   PetscFunctionReturn(0);
983 }
984 
985 #undef __FUNCT__
986 #define __FUNCT__ "TSGetTimeStep"
987 /*@
988    TSGetTimeStep - Gets the current timestep size.
989 
990    Not Collective
991 
992    Input Parameter:
993 .  ts - the TS context obtained from TSCreate()
994 
995    Output Parameter:
996 .  dt - the current timestep size
997 
998    Level: intermediate
999 
1000 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1001 
1002 .keywords: TS, get, timestep
1003 @*/
1004 PetscErrorCode  TSGetTimeStep(TS ts,PetscReal* dt)
1005 {
1006   PetscFunctionBegin;
1007   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1008   PetscValidDoublePointer(dt,2);
1009   *dt = ts->time_step;
1010   PetscFunctionReturn(0);
1011 }
1012 
1013 #undef __FUNCT__
1014 #define __FUNCT__ "TSGetSolution"
1015 /*@
1016    TSGetSolution - Returns the solution at the present timestep. It
1017    is valid to call this routine inside the function that you are evaluating
1018    in order to move to the new timestep. This vector not changed until
1019    the solution at the next timestep has been calculated.
1020 
1021    Not Collective, but Vec returned is parallel if TS is parallel
1022 
1023    Input Parameter:
1024 .  ts - the TS context obtained from TSCreate()
1025 
1026    Output Parameter:
1027 .  v - the vector containing the solution
1028 
1029    Level: intermediate
1030 
1031 .seealso: TSGetTimeStep()
1032 
1033 .keywords: TS, timestep, get, solution
1034 @*/
1035 PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1036 {
1037   PetscFunctionBegin;
1038   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1039   PetscValidPointer(v,2);
1040   *v = ts->vec_sol;
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 /* ----- Routines to initialize and destroy a timestepper ---- */
1045 #undef __FUNCT__
1046 #define __FUNCT__ "TSSetProblemType"
1047 /*@
1048   TSSetProblemType - Sets the type of problem to be solved.
1049 
1050   Not collective
1051 
1052   Input Parameters:
1053 + ts   - The TS
1054 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1055 .vb
1056          U_t = A U
1057          U_t = A(t) U
1058          U_t = F(t,U)
1059 .ve
1060 
1061    Level: beginner
1062 
1063 .keywords: TS, problem type
1064 .seealso: TSSetUp(), TSProblemType, TS
1065 @*/
1066 PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1067 {
1068   PetscErrorCode ierr;
1069 
1070   PetscFunctionBegin;
1071   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1072   ts->problem_type = type;
1073   if (type == TS_LINEAR) {
1074     SNES snes;
1075     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1076     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
1077   }
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 #undef __FUNCT__
1082 #define __FUNCT__ "TSGetProblemType"
1083 /*@C
1084   TSGetProblemType - Gets the type of problem to be solved.
1085 
1086   Not collective
1087 
1088   Input Parameter:
1089 . ts   - The TS
1090 
1091   Output Parameter:
1092 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1093 .vb
1094          M U_t = A U
1095          M(t) U_t = A(t) U
1096          U_t = F(t,U)
1097 .ve
1098 
1099    Level: beginner
1100 
1101 .keywords: TS, problem type
1102 .seealso: TSSetUp(), TSProblemType, TS
1103 @*/
1104 PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1105 {
1106   PetscFunctionBegin;
1107   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1108   PetscValidIntPointer(type,2);
1109   *type = ts->problem_type;
1110   PetscFunctionReturn(0);
1111 }
1112 
1113 #undef __FUNCT__
1114 #define __FUNCT__ "TSSetUp"
1115 /*@
1116    TSSetUp - Sets up the internal data structures for the later use
1117    of a timestepper.
1118 
1119    Collective on TS
1120 
1121    Input Parameter:
1122 .  ts - the TS context obtained from TSCreate()
1123 
1124    Notes:
1125    For basic use of the TS solvers the user need not explicitly call
1126    TSSetUp(), since these actions will automatically occur during
1127    the call to TSStep().  However, if one wishes to control this
1128    phase separately, TSSetUp() should be called after TSCreate()
1129    and optional routines of the form TSSetXXX(), but before TSStep().
1130 
1131    Level: advanced
1132 
1133 .keywords: TS, timestep, setup
1134 
1135 .seealso: TSCreate(), TSStep(), TSDestroy()
1136 @*/
1137 PetscErrorCode  TSSetUp(TS ts)
1138 {
1139   PetscErrorCode ierr;
1140 
1141   PetscFunctionBegin;
1142   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1143   if (ts->setupcalled) PetscFunctionReturn(0);
1144 
1145   if (!((PetscObject)ts)->type_name) {
1146     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1147   }
1148 
1149   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1150 
1151   if (ts->ops->setup) {
1152     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1153   }
1154 
1155   ts->setupcalled = PETSC_TRUE;
1156   PetscFunctionReturn(0);
1157 }
1158 
1159 #undef __FUNCT__
1160 #define __FUNCT__ "TSReset"
1161 /*@
1162    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1163 
1164    Collective on TS
1165 
1166    Input Parameter:
1167 .  ts - the TS context obtained from TSCreate()
1168 
1169    Level: beginner
1170 
1171 .keywords: TS, timestep, reset
1172 
1173 .seealso: TSCreate(), TSSetup(), TSDestroy()
1174 @*/
1175 PetscErrorCode  TSReset(TS ts)
1176 {
1177   PetscErrorCode ierr;
1178 
1179   PetscFunctionBegin;
1180   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1181   if (ts->ops->reset) {
1182     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1183   }
1184   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
1185   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
1186   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1187   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
1188   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1189   if (ts->work) {ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);}
1190   ts->setupcalled = PETSC_FALSE;
1191   PetscFunctionReturn(0);
1192 }
1193 
1194 #undef __FUNCT__
1195 #define __FUNCT__ "TSDestroy"
1196 /*@
1197    TSDestroy - Destroys the timestepper context that was created
1198    with TSCreate().
1199 
1200    Collective on TS
1201 
1202    Input Parameter:
1203 .  ts - the TS context obtained from TSCreate()
1204 
1205    Level: beginner
1206 
1207 .keywords: TS, timestepper, destroy
1208 
1209 .seealso: TSCreate(), TSSetUp(), TSSolve()
1210 @*/
1211 PetscErrorCode  TSDestroy(TS *ts)
1212 {
1213   PetscErrorCode ierr;
1214 
1215   PetscFunctionBegin;
1216   if (!*ts) PetscFunctionReturn(0);
1217   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
1218   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1219 
1220   ierr = TSReset((*ts));CHKERRQ(ierr);
1221 
1222   /* if memory was published with AMS then destroy it */
1223   ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr);
1224   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
1225 
1226   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
1227   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
1228   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
1229 
1230   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1231   PetscFunctionReturn(0);
1232 }
1233 
1234 #undef __FUNCT__
1235 #define __FUNCT__ "TSGetSNES"
1236 /*@
1237    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1238    a TS (timestepper) context. Valid only for nonlinear problems.
1239 
1240    Not Collective, but SNES is parallel if TS is parallel
1241 
1242    Input Parameter:
1243 .  ts - the TS context obtained from TSCreate()
1244 
1245    Output Parameter:
1246 .  snes - the nonlinear solver context
1247 
1248    Notes:
1249    The user can then directly manipulate the SNES context to set various
1250    options, etc.  Likewise, the user can then extract and manipulate the
1251    KSP, KSP, and PC contexts as well.
1252 
1253    TSGetSNES() does not work for integrators that do not use SNES; in
1254    this case TSGetSNES() returns PETSC_NULL in snes.
1255 
1256    Level: beginner
1257 
1258 .keywords: timestep, get, SNES
1259 @*/
1260 PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1261 {
1262   PetscErrorCode ierr;
1263 
1264   PetscFunctionBegin;
1265   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1266   PetscValidPointer(snes,2);
1267   if (!ts->snes) {
1268     ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
1269     ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr);
1270     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
1271     if (ts->problem_type == TS_LINEAR) {
1272       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
1273     }
1274   }
1275   *snes = ts->snes;
1276   PetscFunctionReturn(0);
1277 }
1278 
1279 #undef __FUNCT__
1280 #define __FUNCT__ "TSGetKSP"
1281 /*@
1282    TSGetKSP - Returns the KSP (linear solver) associated with
1283    a TS (timestepper) context.
1284 
1285    Not Collective, but KSP is parallel if TS is parallel
1286 
1287    Input Parameter:
1288 .  ts - the TS context obtained from TSCreate()
1289 
1290    Output Parameter:
1291 .  ksp - the nonlinear solver context
1292 
1293    Notes:
1294    The user can then directly manipulate the KSP context to set various
1295    options, etc.  Likewise, the user can then extract and manipulate the
1296    KSP and PC contexts as well.
1297 
1298    TSGetKSP() does not work for integrators that do not use KSP;
1299    in this case TSGetKSP() returns PETSC_NULL in ksp.
1300 
1301    Level: beginner
1302 
1303 .keywords: timestep, get, KSP
1304 @*/
1305 PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1306 {
1307   PetscErrorCode ierr;
1308   SNES           snes;
1309 
1310   PetscFunctionBegin;
1311   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1312   PetscValidPointer(ksp,2);
1313   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1314   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1315   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1316   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
1317   PetscFunctionReturn(0);
1318 }
1319 
1320 /* ----------- Routines to set solver parameters ---------- */
1321 
1322 #undef __FUNCT__
1323 #define __FUNCT__ "TSGetDuration"
1324 /*@
1325    TSGetDuration - Gets the maximum number of timesteps to use and
1326    maximum time for iteration.
1327 
1328    Not Collective
1329 
1330    Input Parameters:
1331 +  ts       - the TS context obtained from TSCreate()
1332 .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1333 -  maxtime  - final time to iterate to, or PETSC_NULL
1334 
1335    Level: intermediate
1336 
1337 .keywords: TS, timestep, get, maximum, iterations, time
1338 @*/
1339 PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1340 {
1341   PetscFunctionBegin;
1342   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1343   if (maxsteps) {
1344     PetscValidIntPointer(maxsteps,2);
1345     *maxsteps = ts->max_steps;
1346   }
1347   if (maxtime ) {
1348     PetscValidScalarPointer(maxtime,3);
1349     *maxtime  = ts->max_time;
1350   }
1351   PetscFunctionReturn(0);
1352 }
1353 
1354 #undef __FUNCT__
1355 #define __FUNCT__ "TSSetDuration"
1356 /*@
1357    TSSetDuration - Sets the maximum number of timesteps to use and
1358    maximum time for iteration.
1359 
1360    Logically Collective on TS
1361 
1362    Input Parameters:
1363 +  ts - the TS context obtained from TSCreate()
1364 .  maxsteps - maximum number of iterations to use
1365 -  maxtime - final time to iterate to
1366 
1367    Options Database Keys:
1368 .  -ts_max_steps <maxsteps> - Sets maxsteps
1369 .  -ts_max_time <maxtime> - Sets maxtime
1370 
1371    Notes:
1372    The default maximum number of iterations is 5000. Default time is 5.0
1373 
1374    Level: intermediate
1375 
1376 .keywords: TS, timestep, set, maximum, iterations
1377 @*/
1378 PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1379 {
1380   PetscFunctionBegin;
1381   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1382   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1383   PetscValidLogicalCollectiveReal(ts,maxtime,2);
1384   ts->max_steps = maxsteps;
1385   ts->max_time  = maxtime;
1386   PetscFunctionReturn(0);
1387 }
1388 
1389 #undef __FUNCT__
1390 #define __FUNCT__ "TSSetSolution"
1391 /*@
1392    TSSetSolution - Sets the initial solution vector
1393    for use by the TS routines.
1394 
1395    Logically Collective on TS and Vec
1396 
1397    Input Parameters:
1398 +  ts - the TS context obtained from TSCreate()
1399 -  x - the solution vector
1400 
1401    Level: beginner
1402 
1403 .keywords: TS, timestep, set, solution, initial conditions
1404 @*/
1405 PetscErrorCode  TSSetSolution(TS ts,Vec x)
1406 {
1407   PetscErrorCode ierr;
1408 
1409   PetscFunctionBegin;
1410   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1411   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1412   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
1413   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1414   ts->vec_sol = x;
1415   PetscFunctionReturn(0);
1416 }
1417 
1418 #undef __FUNCT__
1419 #define __FUNCT__ "TSSetPreStep"
1420 /*@C
1421   TSSetPreStep - Sets the general-purpose function
1422   called once at the beginning of each time step.
1423 
1424   Logically Collective on TS
1425 
1426   Input Parameters:
1427 + ts   - The TS context obtained from TSCreate()
1428 - func - The function
1429 
1430   Calling sequence of func:
1431 . func (TS ts);
1432 
1433   Level: intermediate
1434 
1435 .keywords: TS, timestep
1436 @*/
1437 PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1438 {
1439   PetscFunctionBegin;
1440   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1441   ts->ops->prestep = func;
1442   PetscFunctionReturn(0);
1443 }
1444 
1445 #undef __FUNCT__
1446 #define __FUNCT__ "TSPreStep"
1447 /*@C
1448   TSPreStep - Runs the user-defined pre-step function.
1449 
1450   Collective on TS
1451 
1452   Input Parameters:
1453 . ts   - The TS context obtained from TSCreate()
1454 
1455   Notes:
1456   TSPreStep() is typically used within time stepping implementations,
1457   so most users would not generally call this routine themselves.
1458 
1459   Level: developer
1460 
1461 .keywords: TS, timestep
1462 @*/
1463 PetscErrorCode  TSPreStep(TS ts)
1464 {
1465   PetscErrorCode ierr;
1466 
1467   PetscFunctionBegin;
1468   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1469   if (ts->ops->prestep) {
1470     PetscStackPush("TS PreStep function");
1471     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1472     PetscStackPop;
1473   }
1474   PetscFunctionReturn(0);
1475 }
1476 
1477 #undef __FUNCT__
1478 #define __FUNCT__ "TSSetPostStep"
1479 /*@C
1480   TSSetPostStep - Sets the general-purpose function
1481   called once at the end of each time step.
1482 
1483   Logically Collective on TS
1484 
1485   Input Parameters:
1486 + ts   - The TS context obtained from TSCreate()
1487 - func - The function
1488 
1489   Calling sequence of func:
1490 . func (TS ts);
1491 
1492   Level: intermediate
1493 
1494 .keywords: TS, timestep
1495 @*/
1496 PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1497 {
1498   PetscFunctionBegin;
1499   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1500   ts->ops->poststep = func;
1501   PetscFunctionReturn(0);
1502 }
1503 
1504 #undef __FUNCT__
1505 #define __FUNCT__ "TSPostStep"
1506 /*@C
1507   TSPostStep - Runs the user-defined post-step function.
1508 
1509   Collective on TS
1510 
1511   Input Parameters:
1512 . ts   - The TS context obtained from TSCreate()
1513 
1514   Notes:
1515   TSPostStep() is typically used within time stepping implementations,
1516   so most users would not generally call this routine themselves.
1517 
1518   Level: developer
1519 
1520 .keywords: TS, timestep
1521 @*/
1522 PetscErrorCode  TSPostStep(TS ts)
1523 {
1524   PetscErrorCode ierr;
1525 
1526   PetscFunctionBegin;
1527   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1528   if (ts->ops->poststep) {
1529     PetscStackPush("TS PostStep function");
1530     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1531     PetscStackPop;
1532   }
1533   PetscFunctionReturn(0);
1534 }
1535 
1536 /* ------------ Routines to set performance monitoring options ----------- */
1537 
1538 #undef __FUNCT__
1539 #define __FUNCT__ "TSMonitorSet"
1540 /*@C
1541    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1542    timestep to display the iteration's  progress.
1543 
1544    Logically Collective on TS
1545 
1546    Input Parameters:
1547 +  ts - the TS context obtained from TSCreate()
1548 .  func - monitoring routine
1549 .  mctx - [optional] user-defined context for private data for the
1550              monitor routine (use PETSC_NULL if no context is desired)
1551 -  monitordestroy - [optional] routine that frees monitor context
1552           (may be PETSC_NULL)
1553 
1554    Calling sequence of func:
1555 $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1556 
1557 +    ts - the TS context
1558 .    steps - iteration number
1559 .    time - current time
1560 .    x - current iterate
1561 -    mctx - [optional] monitoring context
1562 
1563    Notes:
1564    This routine adds an additional monitor to the list of monitors that
1565    already has been loaded.
1566 
1567    Fortran notes: Only a single monitor function can be set for each TS object
1568 
1569    Level: intermediate
1570 
1571 .keywords: TS, timestep, set, monitor
1572 
1573 .seealso: TSMonitorDefault(), TSMonitorCancel()
1574 @*/
1575 PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1576 {
1577   PetscFunctionBegin;
1578   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1579   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1580   ts->monitor[ts->numbermonitors]           = monitor;
1581   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1582   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1583   PetscFunctionReturn(0);
1584 }
1585 
1586 #undef __FUNCT__
1587 #define __FUNCT__ "TSMonitorCancel"
1588 /*@C
1589    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1590 
1591    Logically Collective on TS
1592 
1593    Input Parameters:
1594 .  ts - the TS context obtained from TSCreate()
1595 
1596    Notes:
1597    There is no way to remove a single, specific monitor.
1598 
1599    Level: intermediate
1600 
1601 .keywords: TS, timestep, set, monitor
1602 
1603 .seealso: TSMonitorDefault(), TSMonitorSet()
1604 @*/
1605 PetscErrorCode  TSMonitorCancel(TS ts)
1606 {
1607   PetscErrorCode ierr;
1608   PetscInt       i;
1609 
1610   PetscFunctionBegin;
1611   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1612   for (i=0; i<ts->numbermonitors; i++) {
1613     if (ts->mdestroy[i]) {
1614       ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
1615     }
1616   }
1617   ts->numbermonitors = 0;
1618   PetscFunctionReturn(0);
1619 }
1620 
1621 #undef __FUNCT__
1622 #define __FUNCT__ "TSMonitorDefault"
1623 /*@
1624    TSMonitorDefault - Sets the Default monitor
1625 
1626    Level: intermediate
1627 
1628 .keywords: TS, set, monitor
1629 
1630 .seealso: TSMonitorDefault(), TSMonitorSet()
1631 @*/
1632 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1633 {
1634   PetscErrorCode ierr;
1635   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1636 
1637   PetscFunctionBegin;
1638   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1639   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1640   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1641   PetscFunctionReturn(0);
1642 }
1643 
1644 #undef __FUNCT__
1645 #define __FUNCT__ "TSStep"
1646 /*@
1647    TSStep - Steps the requested number of timesteps.
1648 
1649    Collective on TS
1650 
1651    Input Parameter:
1652 .  ts - the TS context obtained from TSCreate()
1653 
1654    Output Parameters:
1655 +  steps - number of iterations until termination
1656 -  ptime - time until termination
1657 
1658    Level: beginner
1659 
1660 .keywords: TS, timestep, solve
1661 
1662 .seealso: TSCreate(), TSSetUp(), TSDestroy()
1663 @*/
1664 PetscErrorCode  TSStep(TS ts)
1665 {
1666   PetscErrorCode ierr;
1667 
1668   PetscFunctionBegin;
1669   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1670 
1671   ierr = TSSetUp(ts);CHKERRQ(ierr);
1672 
1673   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1674   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
1675   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1676   PetscFunctionReturn(0);
1677 }
1678 
1679 #undef __FUNCT__
1680 #define __FUNCT__ "TSSolve"
1681 /*@
1682    TSSolve - Steps the requested number of timesteps.
1683 
1684    Collective on TS
1685 
1686    Input Parameter:
1687 +  ts - the TS context obtained from TSCreate()
1688 -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
1689 
1690    Level: beginner
1691 
1692 .keywords: TS, timestep, solve
1693 
1694 .seealso: TSCreate(), TSSetSolution(), TSStep()
1695 @*/
1696 PetscErrorCode  TSSolve(TS ts, Vec x)
1697 {
1698   PetscInt       i;
1699   PetscErrorCode ierr;
1700 
1701   PetscFunctionBegin;
1702   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1703   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1704   ierr = TSSetSolution(ts,x); CHKERRQ(ierr);
1705   /* reset time step and iteration counters */
1706   ts->steps = 0;
1707   ts->linear_its = 0;
1708   ts->nonlinear_its = 0;
1709   ts->reason = TS_CONVERGED_ITERATING;
1710   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1711 
1712   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
1713     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
1714   } else {
1715     /* steps the requested number of timesteps. */
1716     for (i=0; !ts->reason; i++) {
1717       ierr = TSPreStep(ts);CHKERRQ(ierr);
1718       ierr = TSStep(ts);CHKERRQ(ierr);
1719       if (ts->reason < 0) {
1720         if (ts->errorifstepfailed) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed");
1721       } else if (i >= ts->max_steps) {
1722         ts->reason = TS_CONVERGED_ITS;
1723       } else if (ts->ptime >= ts->max_time) {
1724         ts->reason = TS_CONVERGED_TIME;
1725       }
1726       ierr = TSPostStep(ts);CHKERRQ(ierr);
1727       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1728     }
1729   }
1730   if (!PetscPreLoadingOn) {
1731     ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr);
1732   }
1733   PetscFunctionReturn(0);
1734 }
1735 
1736 #undef __FUNCT__
1737 #define __FUNCT__ "TSMonitor"
1738 /*
1739      Runs the user provided monitor routines, if they exists.
1740 */
1741 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1742 {
1743   PetscErrorCode ierr;
1744   PetscInt       i,n = ts->numbermonitors;
1745 
1746   PetscFunctionBegin;
1747   for (i=0; i<n; i++) {
1748     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1749   }
1750   PetscFunctionReturn(0);
1751 }
1752 
1753 /* ------------------------------------------------------------------------*/
1754 
1755 #undef __FUNCT__
1756 #define __FUNCT__ "TSMonitorLGCreate"
1757 /*@C
1758    TSMonitorLGCreate - Creates a line graph context for use with
1759    TS to monitor convergence of preconditioned residual norms.
1760 
1761    Collective on TS
1762 
1763    Input Parameters:
1764 +  host - the X display to open, or null for the local machine
1765 .  label - the title to put in the title bar
1766 .  x, y - the screen coordinates of the upper left coordinate of the window
1767 -  m, n - the screen width and height in pixels
1768 
1769    Output Parameter:
1770 .  draw - the drawing context
1771 
1772    Options Database Key:
1773 .  -ts_monitor_draw - automatically sets line graph monitor
1774 
1775    Notes:
1776    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1777 
1778    Level: intermediate
1779 
1780 .keywords: TS, monitor, line graph, residual, seealso
1781 
1782 .seealso: TSMonitorLGDestroy(), TSMonitorSet()
1783 
1784 @*/
1785 PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1786 {
1787   PetscDraw      win;
1788   PetscErrorCode ierr;
1789 
1790   PetscFunctionBegin;
1791   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1792   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1793   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1794   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1795 
1796   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1797   PetscFunctionReturn(0);
1798 }
1799 
1800 #undef __FUNCT__
1801 #define __FUNCT__ "TSMonitorLG"
1802 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1803 {
1804   PetscDrawLG    lg = (PetscDrawLG) monctx;
1805   PetscReal      x,y = ptime;
1806   PetscErrorCode ierr;
1807 
1808   PetscFunctionBegin;
1809   if (!monctx) {
1810     MPI_Comm    comm;
1811     PetscViewer viewer;
1812 
1813     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1814     viewer = PETSC_VIEWER_DRAW_(comm);
1815     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
1816   }
1817 
1818   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1819   x = (PetscReal)n;
1820   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1821   if (n < 20 || (n % 5)) {
1822     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1823   }
1824   PetscFunctionReturn(0);
1825 }
1826 
1827 #undef __FUNCT__
1828 #define __FUNCT__ "TSMonitorLGDestroy"
1829 /*@C
1830    TSMonitorLGDestroy - Destroys a line graph context that was created
1831    with TSMonitorLGCreate().
1832 
1833    Collective on PetscDrawLG
1834 
1835    Input Parameter:
1836 .  draw - the drawing context
1837 
1838    Level: intermediate
1839 
1840 .keywords: TS, monitor, line graph, destroy
1841 
1842 .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1843 @*/
1844 PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
1845 {
1846   PetscDraw      draw;
1847   PetscErrorCode ierr;
1848 
1849   PetscFunctionBegin;
1850   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
1851   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
1852   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1853   PetscFunctionReturn(0);
1854 }
1855 
1856 #undef __FUNCT__
1857 #define __FUNCT__ "TSGetTime"
1858 /*@
1859    TSGetTime - Gets the current time.
1860 
1861    Not Collective
1862 
1863    Input Parameter:
1864 .  ts - the TS context obtained from TSCreate()
1865 
1866    Output Parameter:
1867 .  t  - the current time
1868 
1869    Level: beginner
1870 
1871 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1872 
1873 .keywords: TS, get, time
1874 @*/
1875 PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
1876 {
1877   PetscFunctionBegin;
1878   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1879   PetscValidDoublePointer(t,2);
1880   *t = ts->ptime;
1881   PetscFunctionReturn(0);
1882 }
1883 
1884 #undef __FUNCT__
1885 #define __FUNCT__ "TSSetTime"
1886 /*@
1887    TSSetTime - Allows one to reset the time.
1888 
1889    Logically Collective on TS
1890 
1891    Input Parameters:
1892 +  ts - the TS context obtained from TSCreate()
1893 -  time - the time
1894 
1895    Level: intermediate
1896 
1897 .seealso: TSGetTime(), TSSetDuration()
1898 
1899 .keywords: TS, set, time
1900 @*/
1901 PetscErrorCode  TSSetTime(TS ts, PetscReal t)
1902 {
1903   PetscFunctionBegin;
1904   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1905   PetscValidLogicalCollectiveReal(ts,t,2);
1906   ts->ptime = t;
1907   PetscFunctionReturn(0);
1908 }
1909 
1910 #undef __FUNCT__
1911 #define __FUNCT__ "TSSetOptionsPrefix"
1912 /*@C
1913    TSSetOptionsPrefix - Sets the prefix used for searching for all
1914    TS options in the database.
1915 
1916    Logically Collective on TS
1917 
1918    Input Parameter:
1919 +  ts     - The TS context
1920 -  prefix - The prefix to prepend to all option names
1921 
1922    Notes:
1923    A hyphen (-) must NOT be given at the beginning of the prefix name.
1924    The first character of all runtime options is AUTOMATICALLY the
1925    hyphen.
1926 
1927    Level: advanced
1928 
1929 .keywords: TS, set, options, prefix, database
1930 
1931 .seealso: TSSetFromOptions()
1932 
1933 @*/
1934 PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
1935 {
1936   PetscErrorCode ierr;
1937   SNES           snes;
1938 
1939   PetscFunctionBegin;
1940   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1941   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1942   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1943   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
1944   PetscFunctionReturn(0);
1945 }
1946 
1947 
1948 #undef __FUNCT__
1949 #define __FUNCT__ "TSAppendOptionsPrefix"
1950 /*@C
1951    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1952    TS options in the database.
1953 
1954    Logically Collective on TS
1955 
1956    Input Parameter:
1957 +  ts     - The TS context
1958 -  prefix - The prefix to prepend to all option names
1959 
1960    Notes:
1961    A hyphen (-) must NOT be given at the beginning of the prefix name.
1962    The first character of all runtime options is AUTOMATICALLY the
1963    hyphen.
1964 
1965    Level: advanced
1966 
1967 .keywords: TS, append, options, prefix, database
1968 
1969 .seealso: TSGetOptionsPrefix()
1970 
1971 @*/
1972 PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
1973 {
1974   PetscErrorCode ierr;
1975   SNES           snes;
1976 
1977   PetscFunctionBegin;
1978   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1979   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1980   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1981   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
1982   PetscFunctionReturn(0);
1983 }
1984 
1985 #undef __FUNCT__
1986 #define __FUNCT__ "TSGetOptionsPrefix"
1987 /*@C
1988    TSGetOptionsPrefix - Sets the prefix used for searching for all
1989    TS options in the database.
1990 
1991    Not Collective
1992 
1993    Input Parameter:
1994 .  ts - The TS context
1995 
1996    Output Parameter:
1997 .  prefix - A pointer to the prefix string used
1998 
1999    Notes: On the fortran side, the user should pass in a string 'prifix' of
2000    sufficient length to hold the prefix.
2001 
2002    Level: intermediate
2003 
2004 .keywords: TS, get, options, prefix, database
2005 
2006 .seealso: TSAppendOptionsPrefix()
2007 @*/
2008 PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2009 {
2010   PetscErrorCode ierr;
2011 
2012   PetscFunctionBegin;
2013   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2014   PetscValidPointer(prefix,2);
2015   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2016   PetscFunctionReturn(0);
2017 }
2018 
2019 #undef __FUNCT__
2020 #define __FUNCT__ "TSGetRHSJacobian"
2021 /*@C
2022    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2023 
2024    Not Collective, but parallel objects are returned if TS is parallel
2025 
2026    Input Parameter:
2027 .  ts  - The TS context obtained from TSCreate()
2028 
2029    Output Parameters:
2030 +  J   - The Jacobian J of F, where U_t = F(U,t)
2031 .  M   - The preconditioner matrix, usually the same as J
2032 .  func - Function to compute the Jacobian of the RHS
2033 -  ctx - User-defined context for Jacobian evaluation routine
2034 
2035    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2036 
2037    Level: intermediate
2038 
2039 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2040 
2041 .keywords: TS, timestep, get, matrix, Jacobian
2042 @*/
2043 PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2044 {
2045   PetscErrorCode ierr;
2046   SNES           snes;
2047 
2048   PetscFunctionBegin;
2049   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2050   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2051   if (func) *func = ts->ops->rhsjacobian;
2052   if (ctx) *ctx = ts->jacP;
2053   PetscFunctionReturn(0);
2054 }
2055 
2056 #undef __FUNCT__
2057 #define __FUNCT__ "TSGetIJacobian"
2058 /*@C
2059    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
2060 
2061    Not Collective, but parallel objects are returned if TS is parallel
2062 
2063    Input Parameter:
2064 .  ts  - The TS context obtained from TSCreate()
2065 
2066    Output Parameters:
2067 +  A   - The Jacobian of F(t,U,U_t)
2068 .  B   - The preconditioner matrix, often the same as A
2069 .  f   - The function to compute the matrices
2070 - ctx - User-defined context for Jacobian evaluation routine
2071 
2072    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2073 
2074    Level: advanced
2075 
2076 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2077 
2078 .keywords: TS, timestep, get, matrix, Jacobian
2079 @*/
2080 PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2081 {
2082   PetscErrorCode ierr;
2083   SNES           snes;
2084 
2085   PetscFunctionBegin;
2086   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2087   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2088   if (f) *f = ts->ops->ijacobian;
2089   if (ctx) *ctx = ts->jacP;
2090   PetscFunctionReturn(0);
2091 }
2092 
2093 #undef __FUNCT__
2094 #define __FUNCT__ "TSMonitorSolution"
2095 /*@C
2096    TSMonitorSolution - Monitors progress of the TS solvers by calling
2097    VecView() for the solution at each timestep
2098 
2099    Collective on TS
2100 
2101    Input Parameters:
2102 +  ts - the TS context
2103 .  step - current time-step
2104 .  ptime - current time
2105 -  dummy - either a viewer or PETSC_NULL
2106 
2107    Level: intermediate
2108 
2109 .keywords: TS,  vector, monitor, view
2110 
2111 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2112 @*/
2113 PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2114 {
2115   PetscErrorCode ierr;
2116   PetscViewer    viewer = (PetscViewer) dummy;
2117 
2118   PetscFunctionBegin;
2119   if (!dummy) {
2120     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2121   }
2122   ierr = VecView(x,viewer);CHKERRQ(ierr);
2123   PetscFunctionReturn(0);
2124 }
2125 
2126 
2127 #undef __FUNCT__
2128 #define __FUNCT__ "TSSetDM"
2129 /*@
2130    TSSetDM - Sets the DM that may be used by some preconditioners
2131 
2132    Logically Collective on TS and DM
2133 
2134    Input Parameters:
2135 +  ts - the preconditioner context
2136 -  dm - the dm
2137 
2138    Level: intermediate
2139 
2140 
2141 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2142 @*/
2143 PetscErrorCode  TSSetDM(TS ts,DM dm)
2144 {
2145   PetscErrorCode ierr;
2146   SNES           snes;
2147 
2148   PetscFunctionBegin;
2149   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2150   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
2151   ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
2152   ts->dm = dm;
2153   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2154   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
2155   PetscFunctionReturn(0);
2156 }
2157 
2158 #undef __FUNCT__
2159 #define __FUNCT__ "TSGetDM"
2160 /*@
2161    TSGetDM - Gets the DM that may be used by some preconditioners
2162 
2163    Not Collective
2164 
2165    Input Parameter:
2166 . ts - the preconditioner context
2167 
2168    Output Parameter:
2169 .  dm - the dm
2170 
2171    Level: intermediate
2172 
2173 
2174 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2175 @*/
2176 PetscErrorCode  TSGetDM(TS ts,DM *dm)
2177 {
2178   PetscFunctionBegin;
2179   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2180   *dm = ts->dm;
2181   PetscFunctionReturn(0);
2182 }
2183 
2184 #undef __FUNCT__
2185 #define __FUNCT__ "SNESTSFormFunction"
2186 /*@
2187    SNESTSFormFunction - Function to evaluate nonlinear residual
2188 
2189    Logically Collective on SNES
2190 
2191    Input Parameter:
2192 + snes - nonlinear solver
2193 . X - the current state at which to evaluate the residual
2194 - ctx - user context, must be a TS
2195 
2196    Output Parameter:
2197 . F - the nonlinear residual
2198 
2199    Notes:
2200    This function is not normally called by users and is automatically registered with the SNES used by TS.
2201    It is most frequently passed to MatFDColoringSetFunction().
2202 
2203    Level: advanced
2204 
2205 .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2206 @*/
2207 PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2208 {
2209   TS ts = (TS)ctx;
2210   PetscErrorCode ierr;
2211 
2212   PetscFunctionBegin;
2213   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2214   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2215   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
2216   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
2217   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
2218   PetscFunctionReturn(0);
2219 }
2220 
2221 #undef __FUNCT__
2222 #define __FUNCT__ "SNESTSFormJacobian"
2223 /*@
2224    SNESTSFormJacobian - Function to evaluate the Jacobian
2225 
2226    Collective on SNES
2227 
2228    Input Parameter:
2229 + snes - nonlinear solver
2230 . X - the current state at which to evaluate the residual
2231 - ctx - user context, must be a TS
2232 
2233    Output Parameter:
2234 + A - the Jacobian
2235 . B - the preconditioning matrix (may be the same as A)
2236 - flag - indicates any structure change in the matrix
2237 
2238    Notes:
2239    This function is not normally called by users and is automatically registered with the SNES used by TS.
2240 
2241    Level: developer
2242 
2243 .seealso: SNESSetJacobian()
2244 @*/
2245 PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2246 {
2247   TS ts = (TS)ctx;
2248   PetscErrorCode ierr;
2249 
2250   PetscFunctionBegin;
2251   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2252   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2253   PetscValidPointer(A,3);
2254   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
2255   PetscValidPointer(B,4);
2256   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
2257   PetscValidPointer(flag,5);
2258   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
2259   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
2260   PetscFunctionReturn(0);
2261 }
2262 
2263 #if defined(PETSC_HAVE_MATLAB_ENGINE)
2264 #include <mex.h>
2265 
2266 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
2267 
2268 #undef __FUNCT__
2269 #define __FUNCT__ "TSComputeFunction_Matlab"
2270 /*
2271    TSComputeFunction_Matlab - Calls the function that has been set with
2272                          TSSetFunctionMatlab().
2273 
2274    Collective on TS
2275 
2276    Input Parameters:
2277 +  snes - the TS context
2278 -  x - input vector
2279 
2280    Output Parameter:
2281 .  y - function vector, as set by TSSetFunction()
2282 
2283    Notes:
2284    TSComputeFunction() is typically used within nonlinear solvers
2285    implementations, so most users would not generally call this routine
2286    themselves.
2287 
2288    Level: developer
2289 
2290 .keywords: TS, nonlinear, compute, function
2291 
2292 .seealso: TSSetFunction(), TSGetFunction()
2293 */
2294 PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
2295 {
2296   PetscErrorCode   ierr;
2297   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2298   int              nlhs = 1,nrhs = 7;
2299   mxArray	   *plhs[1],*prhs[7];
2300   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
2301 
2302   PetscFunctionBegin;
2303   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2304   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2305   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
2306   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
2307   PetscCheckSameComm(snes,1,x,3);
2308   PetscCheckSameComm(snes,1,y,5);
2309 
2310   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2311   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2312   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
2313   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
2314   prhs[0] =  mxCreateDoubleScalar((double)ls);
2315   prhs[1] =  mxCreateDoubleScalar(time);
2316   prhs[2] =  mxCreateDoubleScalar((double)lx);
2317   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2318   prhs[4] =  mxCreateDoubleScalar((double)ly);
2319   prhs[5] =  mxCreateString(sctx->funcname);
2320   prhs[6] =  sctx->ctx;
2321   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
2322   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2323   mxDestroyArray(prhs[0]);
2324   mxDestroyArray(prhs[1]);
2325   mxDestroyArray(prhs[2]);
2326   mxDestroyArray(prhs[3]);
2327   mxDestroyArray(prhs[4]);
2328   mxDestroyArray(prhs[5]);
2329   mxDestroyArray(plhs[0]);
2330   PetscFunctionReturn(0);
2331 }
2332 
2333 
2334 #undef __FUNCT__
2335 #define __FUNCT__ "TSSetFunctionMatlab"
2336 /*
2337    TSSetFunctionMatlab - Sets the function evaluation routine and function
2338    vector for use by the TS routines in solving ODEs
2339    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
2340 
2341    Logically Collective on TS
2342 
2343    Input Parameters:
2344 +  ts - the TS context
2345 -  func - function evaluation routine
2346 
2347    Calling sequence of func:
2348 $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
2349 
2350    Level: beginner
2351 
2352 .keywords: TS, nonlinear, set, function
2353 
2354 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2355 */
2356 PetscErrorCode  TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx)
2357 {
2358   PetscErrorCode  ierr;
2359   TSMatlabContext *sctx;
2360 
2361   PetscFunctionBegin;
2362   /* currently sctx is memory bleed */
2363   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2364   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2365   /*
2366      This should work, but it doesn't
2367   sctx->ctx = ctx;
2368   mexMakeArrayPersistent(sctx->ctx);
2369   */
2370   sctx->ctx = mxDuplicateArray(ctx);
2371   ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
2372   PetscFunctionReturn(0);
2373 }
2374 
2375 #undef __FUNCT__
2376 #define __FUNCT__ "TSComputeJacobian_Matlab"
2377 /*
2378    TSComputeJacobian_Matlab - Calls the function that has been set with
2379                          TSSetJacobianMatlab().
2380 
2381    Collective on TS
2382 
2383    Input Parameters:
2384 +  snes - the TS context
2385 .  x - input vector
2386 .  A, B - the matrices
2387 -  ctx - user context
2388 
2389    Output Parameter:
2390 .  flag - structure of the matrix
2391 
2392    Level: developer
2393 
2394 .keywords: TS, nonlinear, compute, function
2395 
2396 .seealso: TSSetFunction(), TSGetFunction()
2397 @*/
2398 PetscErrorCode  TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
2399 {
2400   PetscErrorCode  ierr;
2401   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2402   int             nlhs = 2,nrhs = 9;
2403   mxArray	  *plhs[2],*prhs[9];
2404   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
2405 
2406   PetscFunctionBegin;
2407   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2408   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2409 
2410   /* call Matlab function in ctx with arguments x and y */
2411 
2412   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2413   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2414   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
2415   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
2416   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
2417   prhs[0] =  mxCreateDoubleScalar((double)ls);
2418   prhs[1] =  mxCreateDoubleScalar((double)time);
2419   prhs[2] =  mxCreateDoubleScalar((double)lx);
2420   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2421   prhs[4] =  mxCreateDoubleScalar((double)shift);
2422   prhs[5] =  mxCreateDoubleScalar((double)lA);
2423   prhs[6] =  mxCreateDoubleScalar((double)lB);
2424   prhs[7] =  mxCreateString(sctx->funcname);
2425   prhs[8] =  sctx->ctx;
2426   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
2427   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2428   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
2429   mxDestroyArray(prhs[0]);
2430   mxDestroyArray(prhs[1]);
2431   mxDestroyArray(prhs[2]);
2432   mxDestroyArray(prhs[3]);
2433   mxDestroyArray(prhs[4]);
2434   mxDestroyArray(prhs[5]);
2435   mxDestroyArray(prhs[6]);
2436   mxDestroyArray(prhs[7]);
2437   mxDestroyArray(plhs[0]);
2438   mxDestroyArray(plhs[1]);
2439   PetscFunctionReturn(0);
2440 }
2441 
2442 
2443 #undef __FUNCT__
2444 #define __FUNCT__ "TSSetJacobianMatlab"
2445 /*
2446    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
2447    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
2448 
2449    Logically Collective on TS
2450 
2451    Input Parameters:
2452 +  snes - the TS context
2453 .  A,B - Jacobian matrices
2454 .  func - function evaluation routine
2455 -  ctx - user context
2456 
2457    Calling sequence of func:
2458 $    flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
2459 
2460 
2461    Level: developer
2462 
2463 .keywords: TS, nonlinear, set, function
2464 
2465 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2466 */
2467 PetscErrorCode  TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx)
2468 {
2469   PetscErrorCode    ierr;
2470   TSMatlabContext *sctx;
2471 
2472   PetscFunctionBegin;
2473   /* currently sctx is memory bleed */
2474   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2475   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2476   /*
2477      This should work, but it doesn't
2478   sctx->ctx = ctx;
2479   mexMakeArrayPersistent(sctx->ctx);
2480   */
2481   sctx->ctx = mxDuplicateArray(ctx);
2482   ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
2483   PetscFunctionReturn(0);
2484 }
2485 
2486 #undef __FUNCT__
2487 #define __FUNCT__ "TSMonitor_Matlab"
2488 /*
2489    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
2490 
2491    Collective on TS
2492 
2493 .seealso: TSSetFunction(), TSGetFunction()
2494 @*/
2495 PetscErrorCode  TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx)
2496 {
2497   PetscErrorCode  ierr;
2498   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2499   int             nlhs = 1,nrhs = 6;
2500   mxArray	  *plhs[1],*prhs[6];
2501   long long int   lx = 0,ls = 0;
2502 
2503   PetscFunctionBegin;
2504   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2505   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
2506 
2507   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2508   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2509   prhs[0] =  mxCreateDoubleScalar((double)ls);
2510   prhs[1] =  mxCreateDoubleScalar((double)it);
2511   prhs[2] =  mxCreateDoubleScalar((double)time);
2512   prhs[3] =  mxCreateDoubleScalar((double)lx);
2513   prhs[4] =  mxCreateString(sctx->funcname);
2514   prhs[5] =  sctx->ctx;
2515   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
2516   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2517   mxDestroyArray(prhs[0]);
2518   mxDestroyArray(prhs[1]);
2519   mxDestroyArray(prhs[2]);
2520   mxDestroyArray(prhs[3]);
2521   mxDestroyArray(prhs[4]);
2522   mxDestroyArray(plhs[0]);
2523   PetscFunctionReturn(0);
2524 }
2525 
2526 
2527 #undef __FUNCT__
2528 #define __FUNCT__ "TSMonitorSetMatlab"
2529 /*
2530    TSMonitorSetMatlab - Sets the monitor function from Matlab
2531 
2532    Level: developer
2533 
2534 .keywords: TS, nonlinear, set, function
2535 
2536 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2537 */
2538 PetscErrorCode  TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx)
2539 {
2540   PetscErrorCode    ierr;
2541   TSMatlabContext *sctx;
2542 
2543   PetscFunctionBegin;
2544   /* currently sctx is memory bleed */
2545   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2546   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2547   /*
2548      This should work, but it doesn't
2549   sctx->ctx = ctx;
2550   mexMakeArrayPersistent(sctx->ctx);
2551   */
2552   sctx->ctx = mxDuplicateArray(ctx);
2553   ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
2554   PetscFunctionReturn(0);
2555 }
2556 
2557 #endif
2558