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