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