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