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