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