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