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