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