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