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