xref: /petsc/src/ts/interface/ts.c (revision b2cd27e8157cb66dc382ecf63346fa09165ec73c)
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   PetscFunctionBegin;
1163   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1164   ts->problem_type = type;
1165   PetscFunctionReturn(0);
1166 }
1167 
1168 #undef __FUNCT__
1169 #define __FUNCT__ "TSGetProblemType"
1170 /*@C
1171   TSGetProblemType - Gets the type of problem to be solved.
1172 
1173   Not collective
1174 
1175   Input Parameter:
1176 . ts   - The TS
1177 
1178   Output Parameter:
1179 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1180 .vb
1181          M U_t = A U
1182          M(t) U_t = A(t) U
1183          U_t = F(t,U)
1184 .ve
1185 
1186    Level: beginner
1187 
1188 .keywords: TS, problem type
1189 .seealso: TSSetUp(), TSProblemType, TS
1190 @*/
1191 PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1192 {
1193   PetscFunctionBegin;
1194   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1195   PetscValidIntPointer(type,2);
1196   *type = ts->problem_type;
1197   PetscFunctionReturn(0);
1198 }
1199 
1200 #undef __FUNCT__
1201 #define __FUNCT__ "TSSetUp"
1202 /*@
1203    TSSetUp - Sets up the internal data structures for the later use
1204    of a timestepper.
1205 
1206    Collective on TS
1207 
1208    Input Parameter:
1209 .  ts - the TS context obtained from TSCreate()
1210 
1211    Notes:
1212    For basic use of the TS solvers the user need not explicitly call
1213    TSSetUp(), since these actions will automatically occur during
1214    the call to TSStep().  However, if one wishes to control this
1215    phase separately, TSSetUp() should be called after TSCreate()
1216    and optional routines of the form TSSetXXX(), but before TSStep().
1217 
1218    Level: advanced
1219 
1220 .keywords: TS, timestep, setup
1221 
1222 .seealso: TSCreate(), TSStep(), TSDestroy()
1223 @*/
1224 PetscErrorCode  TSSetUp(TS ts)
1225 {
1226   PetscErrorCode ierr;
1227 
1228   PetscFunctionBegin;
1229   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1230   if (ts->setupcalled) PetscFunctionReturn(0);
1231 
1232   if (!((PetscObject)ts)->type_name) {
1233     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1234   }
1235 
1236   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1237 
1238   if (ts->ops->setup) {
1239     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1240   }
1241 
1242   ts->setupcalled = PETSC_TRUE;
1243   PetscFunctionReturn(0);
1244 }
1245 
1246 #undef __FUNCT__
1247 #define __FUNCT__ "TSReset"
1248 /*@
1249    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1250 
1251    Collective on TS
1252 
1253    Input Parameter:
1254 .  ts - the TS context obtained from TSCreate()
1255 
1256    Level: beginner
1257 
1258 .keywords: TS, timestep, reset
1259 
1260 .seealso: TSCreate(), TSSetup(), TSDestroy()
1261 @*/
1262 PetscErrorCode  TSReset(TS ts)
1263 {
1264   PetscErrorCode ierr;
1265 
1266   PetscFunctionBegin;
1267   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1268   if (ts->ops->reset) {
1269     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1270   }
1271   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
1272   ierr = MatDestroy(&ts->Alhs);CHKERRQ(ierr);
1273   ierr = MatDestroy(&ts->Blhs);CHKERRQ(ierr);
1274   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1275   if (ts->work) {ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);}
1276   ts->setupcalled = PETSC_FALSE;
1277   PetscFunctionReturn(0);
1278 }
1279 
1280 #undef __FUNCT__
1281 #define __FUNCT__ "TSDestroy"
1282 /*@
1283    TSDestroy - Destroys the timestepper context that was created
1284    with TSCreate().
1285 
1286    Collective on TS
1287 
1288    Input Parameter:
1289 .  ts - the TS context obtained from TSCreate()
1290 
1291    Level: beginner
1292 
1293 .keywords: TS, timestepper, destroy
1294 
1295 .seealso: TSCreate(), TSSetUp(), TSSolve()
1296 @*/
1297 PetscErrorCode  TSDestroy(TS *ts)
1298 {
1299   PetscErrorCode ierr;
1300 
1301   PetscFunctionBegin;
1302   if (!*ts) PetscFunctionReturn(0);
1303   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
1304   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1305 
1306   ierr = TSReset((*ts));CHKERRQ(ierr);
1307 
1308   /* if memory was published with AMS then destroy it */
1309   ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr);
1310   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
1311 
1312   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
1313   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
1314   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
1315 
1316   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1317   PetscFunctionReturn(0);
1318 }
1319 
1320 #undef __FUNCT__
1321 #define __FUNCT__ "TSGetSNES"
1322 /*@
1323    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1324    a TS (timestepper) context. Valid only for nonlinear problems.
1325 
1326    Not Collective, but SNES is parallel if TS is parallel
1327 
1328    Input Parameter:
1329 .  ts - the TS context obtained from TSCreate()
1330 
1331    Output Parameter:
1332 .  snes - the nonlinear solver context
1333 
1334    Notes:
1335    The user can then directly manipulate the SNES context to set various
1336    options, etc.  Likewise, the user can then extract and manipulate the
1337    KSP, KSP, and PC contexts as well.
1338 
1339    TSGetSNES() does not work for integrators that do not use SNES; in
1340    this case TSGetSNES() returns PETSC_NULL in snes.
1341 
1342    Level: beginner
1343 
1344 .keywords: timestep, get, SNES
1345 @*/
1346 PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1347 {
1348   PetscErrorCode ierr;
1349 
1350   PetscFunctionBegin;
1351   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1352   PetscValidPointer(snes,2);
1353   if (!ts->snes) {
1354     ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
1355     ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr);
1356     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
1357   }
1358   *snes = ts->snes;
1359   PetscFunctionReturn(0);
1360 }
1361 
1362 #undef __FUNCT__
1363 #define __FUNCT__ "TSGetKSP"
1364 /*@
1365    TSGetKSP - Returns the KSP (linear solver) associated with
1366    a TS (timestepper) context.
1367 
1368    Not Collective, but KSP is parallel if TS is parallel
1369 
1370    Input Parameter:
1371 .  ts - the TS context obtained from TSCreate()
1372 
1373    Output Parameter:
1374 .  ksp - the nonlinear solver context
1375 
1376    Notes:
1377    The user can then directly manipulate the KSP context to set various
1378    options, etc.  Likewise, the user can then extract and manipulate the
1379    KSP and PC contexts as well.
1380 
1381    TSGetKSP() does not work for integrators that do not use KSP;
1382    in this case TSGetKSP() returns PETSC_NULL in ksp.
1383 
1384    Level: beginner
1385 
1386 .keywords: timestep, get, KSP
1387 @*/
1388 PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1389 {
1390   PetscErrorCode ierr;
1391   SNES           snes;
1392 
1393   PetscFunctionBegin;
1394   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1395   PetscValidPointer(ksp,2);
1396   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1397   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1398   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1399   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
1400   PetscFunctionReturn(0);
1401 }
1402 
1403 /* ----------- Routines to set solver parameters ---------- */
1404 
1405 #undef __FUNCT__
1406 #define __FUNCT__ "TSGetDuration"
1407 /*@
1408    TSGetDuration - Gets the maximum number of timesteps to use and
1409    maximum time for iteration.
1410 
1411    Not Collective
1412 
1413    Input Parameters:
1414 +  ts       - the TS context obtained from TSCreate()
1415 .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1416 -  maxtime  - final time to iterate to, or PETSC_NULL
1417 
1418    Level: intermediate
1419 
1420 .keywords: TS, timestep, get, maximum, iterations, time
1421 @*/
1422 PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1423 {
1424   PetscFunctionBegin;
1425   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1426   if (maxsteps) {
1427     PetscValidIntPointer(maxsteps,2);
1428     *maxsteps = ts->max_steps;
1429   }
1430   if (maxtime ) {
1431     PetscValidScalarPointer(maxtime,3);
1432     *maxtime  = ts->max_time;
1433   }
1434   PetscFunctionReturn(0);
1435 }
1436 
1437 #undef __FUNCT__
1438 #define __FUNCT__ "TSSetDuration"
1439 /*@
1440    TSSetDuration - Sets the maximum number of timesteps to use and
1441    maximum time for iteration.
1442 
1443    Logically Collective on TS
1444 
1445    Input Parameters:
1446 +  ts - the TS context obtained from TSCreate()
1447 .  maxsteps - maximum number of iterations to use
1448 -  maxtime - final time to iterate to
1449 
1450    Options Database Keys:
1451 .  -ts_max_steps <maxsteps> - Sets maxsteps
1452 .  -ts_max_time <maxtime> - Sets maxtime
1453 
1454    Notes:
1455    The default maximum number of iterations is 5000. Default time is 5.0
1456 
1457    Level: intermediate
1458 
1459 .keywords: TS, timestep, set, maximum, iterations
1460 @*/
1461 PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1462 {
1463   PetscFunctionBegin;
1464   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1465   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1466   PetscValidLogicalCollectiveReal(ts,maxtime,2);
1467   ts->max_steps = maxsteps;
1468   ts->max_time  = maxtime;
1469   PetscFunctionReturn(0);
1470 }
1471 
1472 #undef __FUNCT__
1473 #define __FUNCT__ "TSSetSolution"
1474 /*@
1475    TSSetSolution - Sets the initial solution vector
1476    for use by the TS routines.
1477 
1478    Logically Collective on TS and Vec
1479 
1480    Input Parameters:
1481 +  ts - the TS context obtained from TSCreate()
1482 -  x - the solution vector
1483 
1484    Level: beginner
1485 
1486 .keywords: TS, timestep, set, solution, initial conditions
1487 @*/
1488 PetscErrorCode  TSSetSolution(TS ts,Vec x)
1489 {
1490   PetscErrorCode ierr;
1491 
1492   PetscFunctionBegin;
1493   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1494   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1495   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
1496   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1497   ts->vec_sol = x;
1498   PetscFunctionReturn(0);
1499 }
1500 
1501 #undef __FUNCT__
1502 #define __FUNCT__ "TSSetPreStep"
1503 /*@C
1504   TSSetPreStep - Sets the general-purpose function
1505   called once at the beginning of each time step.
1506 
1507   Logically Collective on TS
1508 
1509   Input Parameters:
1510 + ts   - The TS context obtained from TSCreate()
1511 - func - The function
1512 
1513   Calling sequence of func:
1514 . func (TS ts);
1515 
1516   Level: intermediate
1517 
1518 .keywords: TS, timestep
1519 @*/
1520 PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1521 {
1522   PetscFunctionBegin;
1523   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1524   ts->ops->prestep = func;
1525   PetscFunctionReturn(0);
1526 }
1527 
1528 #undef __FUNCT__
1529 #define __FUNCT__ "TSPreStep"
1530 /*@C
1531   TSPreStep - Runs the user-defined pre-step function.
1532 
1533   Collective on TS
1534 
1535   Input Parameters:
1536 . ts   - The TS context obtained from TSCreate()
1537 
1538   Notes:
1539   TSPreStep() is typically used within time stepping implementations,
1540   so most users would not generally call this routine themselves.
1541 
1542   Level: developer
1543 
1544 .keywords: TS, timestep
1545 @*/
1546 PetscErrorCode  TSPreStep(TS ts)
1547 {
1548   PetscErrorCode ierr;
1549 
1550   PetscFunctionBegin;
1551   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1552   if (ts->ops->prestep) {
1553     PetscStackPush("TS PreStep function");
1554     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1555     PetscStackPop;
1556   }
1557   PetscFunctionReturn(0);
1558 }
1559 
1560 #undef __FUNCT__
1561 #define __FUNCT__ "TSSetPostStep"
1562 /*@C
1563   TSSetPostStep - Sets the general-purpose function
1564   called once at the end of each time step.
1565 
1566   Logically Collective on TS
1567 
1568   Input Parameters:
1569 + ts   - The TS context obtained from TSCreate()
1570 - func - The function
1571 
1572   Calling sequence of func:
1573 . func (TS ts);
1574 
1575   Level: intermediate
1576 
1577 .keywords: TS, timestep
1578 @*/
1579 PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1580 {
1581   PetscFunctionBegin;
1582   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1583   ts->ops->poststep = func;
1584   PetscFunctionReturn(0);
1585 }
1586 
1587 #undef __FUNCT__
1588 #define __FUNCT__ "TSPostStep"
1589 /*@C
1590   TSPostStep - Runs the user-defined post-step function.
1591 
1592   Collective on TS
1593 
1594   Input Parameters:
1595 . ts   - The TS context obtained from TSCreate()
1596 
1597   Notes:
1598   TSPostStep() is typically used within time stepping implementations,
1599   so most users would not generally call this routine themselves.
1600 
1601   Level: developer
1602 
1603 .keywords: TS, timestep
1604 @*/
1605 PetscErrorCode  TSPostStep(TS ts)
1606 {
1607   PetscErrorCode ierr;
1608 
1609   PetscFunctionBegin;
1610   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1611   if (ts->ops->poststep) {
1612     PetscStackPush("TS PostStep function");
1613     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1614     PetscStackPop;
1615   }
1616   PetscFunctionReturn(0);
1617 }
1618 
1619 /* ------------ Routines to set performance monitoring options ----------- */
1620 
1621 #undef __FUNCT__
1622 #define __FUNCT__ "TSMonitorSet"
1623 /*@C
1624    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1625    timestep to display the iteration's  progress.
1626 
1627    Logically Collective on TS
1628 
1629    Input Parameters:
1630 +  ts - the TS context obtained from TSCreate()
1631 .  func - monitoring routine
1632 .  mctx - [optional] user-defined context for private data for the
1633              monitor routine (use PETSC_NULL if no context is desired)
1634 -  monitordestroy - [optional] routine that frees monitor context
1635           (may be PETSC_NULL)
1636 
1637    Calling sequence of func:
1638 $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1639 
1640 +    ts - the TS context
1641 .    steps - iteration number
1642 .    time - current time
1643 .    x - current iterate
1644 -    mctx - [optional] monitoring context
1645 
1646    Notes:
1647    This routine adds an additional monitor to the list of monitors that
1648    already has been loaded.
1649 
1650    Fortran notes: Only a single monitor function can be set for each TS object
1651 
1652    Level: intermediate
1653 
1654 .keywords: TS, timestep, set, monitor
1655 
1656 .seealso: TSMonitorDefault(), TSMonitorCancel()
1657 @*/
1658 PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1659 {
1660   PetscFunctionBegin;
1661   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1662   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1663   ts->monitor[ts->numbermonitors]           = monitor;
1664   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1665   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1666   PetscFunctionReturn(0);
1667 }
1668 
1669 #undef __FUNCT__
1670 #define __FUNCT__ "TSMonitorCancel"
1671 /*@C
1672    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1673 
1674    Logically Collective on TS
1675 
1676    Input Parameters:
1677 .  ts - the TS context obtained from TSCreate()
1678 
1679    Notes:
1680    There is no way to remove a single, specific monitor.
1681 
1682    Level: intermediate
1683 
1684 .keywords: TS, timestep, set, monitor
1685 
1686 .seealso: TSMonitorDefault(), TSMonitorSet()
1687 @*/
1688 PetscErrorCode  TSMonitorCancel(TS ts)
1689 {
1690   PetscErrorCode ierr;
1691   PetscInt       i;
1692 
1693   PetscFunctionBegin;
1694   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1695   for (i=0; i<ts->numbermonitors; i++) {
1696     if (ts->mdestroy[i]) {
1697       ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
1698     }
1699   }
1700   ts->numbermonitors = 0;
1701   PetscFunctionReturn(0);
1702 }
1703 
1704 #undef __FUNCT__
1705 #define __FUNCT__ "TSMonitorDefault"
1706 /*@
1707    TSMonitorDefault - Sets the Default monitor
1708 
1709    Level: intermediate
1710 
1711 .keywords: TS, set, monitor
1712 
1713 .seealso: TSMonitorDefault(), TSMonitorSet()
1714 @*/
1715 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1716 {
1717   PetscErrorCode ierr;
1718   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1719 
1720   PetscFunctionBegin;
1721   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1722   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1723   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1724   PetscFunctionReturn(0);
1725 }
1726 
1727 #undef __FUNCT__
1728 #define __FUNCT__ "TSStep"
1729 /*@
1730    TSStep - Steps the requested number of timesteps.
1731 
1732    Collective on TS
1733 
1734    Input Parameter:
1735 .  ts - the TS context obtained from TSCreate()
1736 
1737    Output Parameters:
1738 +  steps - number of iterations until termination
1739 -  ptime - time until termination
1740 
1741    Level: beginner
1742 
1743 .keywords: TS, timestep, solve
1744 
1745 .seealso: TSCreate(), TSSetUp(), TSDestroy()
1746 @*/
1747 PetscErrorCode  TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1748 {
1749   PetscErrorCode ierr;
1750 
1751   PetscFunctionBegin;
1752   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1753 
1754   ierr = TSSetUp(ts);CHKERRQ(ierr);
1755 
1756   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1757   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1758   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1759 
1760   if (!PetscPreLoadingOn) {
1761     ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr);
1762   }
1763   PetscFunctionReturn(0);
1764 }
1765 
1766 #undef __FUNCT__
1767 #define __FUNCT__ "TSSolve"
1768 /*@
1769    TSSolve - Steps the requested number of timesteps.
1770 
1771    Collective on TS
1772 
1773    Input Parameter:
1774 +  ts - the TS context obtained from TSCreate()
1775 -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
1776 
1777    Level: beginner
1778 
1779 .keywords: TS, timestep, solve
1780 
1781 .seealso: TSCreate(), TSSetSolution(), TSStep()
1782 @*/
1783 PetscErrorCode  TSSolve(TS ts, Vec x)
1784 {
1785   PetscInt       steps;
1786   PetscReal      ptime;
1787   PetscErrorCode ierr;
1788 
1789   PetscFunctionBegin;
1790   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1791   /* set solution vector if provided */
1792   if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); }
1793   /* reset time step and iteration counters */
1794   ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0;
1795   /* steps the requested number of timesteps. */
1796   ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr);
1797   PetscFunctionReturn(0);
1798 }
1799 
1800 #undef __FUNCT__
1801 #define __FUNCT__ "TSMonitor"
1802 /*
1803      Runs the user provided monitor routines, if they exists.
1804 */
1805 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1806 {
1807   PetscErrorCode ierr;
1808   PetscInt       i,n = ts->numbermonitors;
1809 
1810   PetscFunctionBegin;
1811   for (i=0; i<n; i++) {
1812     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1813   }
1814   PetscFunctionReturn(0);
1815 }
1816 
1817 /* ------------------------------------------------------------------------*/
1818 
1819 #undef __FUNCT__
1820 #define __FUNCT__ "TSMonitorLGCreate"
1821 /*@C
1822    TSMonitorLGCreate - Creates a line graph context for use with
1823    TS to monitor convergence of preconditioned residual norms.
1824 
1825    Collective on TS
1826 
1827    Input Parameters:
1828 +  host - the X display to open, or null for the local machine
1829 .  label - the title to put in the title bar
1830 .  x, y - the screen coordinates of the upper left coordinate of the window
1831 -  m, n - the screen width and height in pixels
1832 
1833    Output Parameter:
1834 .  draw - the drawing context
1835 
1836    Options Database Key:
1837 .  -ts_monitor_draw - automatically sets line graph monitor
1838 
1839    Notes:
1840    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1841 
1842    Level: intermediate
1843 
1844 .keywords: TS, monitor, line graph, residual, seealso
1845 
1846 .seealso: TSMonitorLGDestroy(), TSMonitorSet()
1847 
1848 @*/
1849 PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1850 {
1851   PetscDraw      win;
1852   PetscErrorCode ierr;
1853 
1854   PetscFunctionBegin;
1855   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1856   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1857   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1858   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1859 
1860   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1861   PetscFunctionReturn(0);
1862 }
1863 
1864 #undef __FUNCT__
1865 #define __FUNCT__ "TSMonitorLG"
1866 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1867 {
1868   PetscDrawLG    lg = (PetscDrawLG) monctx;
1869   PetscReal      x,y = ptime;
1870   PetscErrorCode ierr;
1871 
1872   PetscFunctionBegin;
1873   if (!monctx) {
1874     MPI_Comm    comm;
1875     PetscViewer viewer;
1876 
1877     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1878     viewer = PETSC_VIEWER_DRAW_(comm);
1879     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
1880   }
1881 
1882   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1883   x = (PetscReal)n;
1884   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1885   if (n < 20 || (n % 5)) {
1886     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1887   }
1888   PetscFunctionReturn(0);
1889 }
1890 
1891 #undef __FUNCT__
1892 #define __FUNCT__ "TSMonitorLGDestroy"
1893 /*@C
1894    TSMonitorLGDestroy - Destroys a line graph context that was created
1895    with TSMonitorLGCreate().
1896 
1897    Collective on PetscDrawLG
1898 
1899    Input Parameter:
1900 .  draw - the drawing context
1901 
1902    Level: intermediate
1903 
1904 .keywords: TS, monitor, line graph, destroy
1905 
1906 .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1907 @*/
1908 PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
1909 {
1910   PetscDraw      draw;
1911   PetscErrorCode ierr;
1912 
1913   PetscFunctionBegin;
1914   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
1915   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
1916   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1917   PetscFunctionReturn(0);
1918 }
1919 
1920 #undef __FUNCT__
1921 #define __FUNCT__ "TSGetTime"
1922 /*@
1923    TSGetTime - Gets the current time.
1924 
1925    Not Collective
1926 
1927    Input Parameter:
1928 .  ts - the TS context obtained from TSCreate()
1929 
1930    Output Parameter:
1931 .  t  - the current time
1932 
1933    Level: beginner
1934 
1935 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1936 
1937 .keywords: TS, get, time
1938 @*/
1939 PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
1940 {
1941   PetscFunctionBegin;
1942   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1943   PetscValidDoublePointer(t,2);
1944   *t = ts->ptime;
1945   PetscFunctionReturn(0);
1946 }
1947 
1948 #undef __FUNCT__
1949 #define __FUNCT__ "TSSetTime"
1950 /*@
1951    TSSetTime - Allows one to reset the time.
1952 
1953    Logically Collective on TS
1954 
1955    Input Parameters:
1956 +  ts - the TS context obtained from TSCreate()
1957 -  time - the time
1958 
1959    Level: intermediate
1960 
1961 .seealso: TSGetTime(), TSSetDuration()
1962 
1963 .keywords: TS, set, time
1964 @*/
1965 PetscErrorCode  TSSetTime(TS ts, PetscReal t)
1966 {
1967   PetscFunctionBegin;
1968   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1969   PetscValidLogicalCollectiveReal(ts,t,2);
1970   ts->ptime = t;
1971   PetscFunctionReturn(0);
1972 }
1973 
1974 #undef __FUNCT__
1975 #define __FUNCT__ "TSSetOptionsPrefix"
1976 /*@C
1977    TSSetOptionsPrefix - Sets the prefix used for searching for all
1978    TS options in the database.
1979 
1980    Logically Collective on TS
1981 
1982    Input Parameter:
1983 +  ts     - The TS context
1984 -  prefix - The prefix to prepend to all option names
1985 
1986    Notes:
1987    A hyphen (-) must NOT be given at the beginning of the prefix name.
1988    The first character of all runtime options is AUTOMATICALLY the
1989    hyphen.
1990 
1991    Level: advanced
1992 
1993 .keywords: TS, set, options, prefix, database
1994 
1995 .seealso: TSSetFromOptions()
1996 
1997 @*/
1998 PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
1999 {
2000   PetscErrorCode ierr;
2001   SNES           snes;
2002 
2003   PetscFunctionBegin;
2004   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2005   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2006   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2007   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2008   PetscFunctionReturn(0);
2009 }
2010 
2011 
2012 #undef __FUNCT__
2013 #define __FUNCT__ "TSAppendOptionsPrefix"
2014 /*@C
2015    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2016    TS options in the database.
2017 
2018    Logically Collective on TS
2019 
2020    Input Parameter:
2021 +  ts     - The TS context
2022 -  prefix - The prefix to prepend to all option names
2023 
2024    Notes:
2025    A hyphen (-) must NOT be given at the beginning of the prefix name.
2026    The first character of all runtime options is AUTOMATICALLY the
2027    hyphen.
2028 
2029    Level: advanced
2030 
2031 .keywords: TS, append, options, prefix, database
2032 
2033 .seealso: TSGetOptionsPrefix()
2034 
2035 @*/
2036 PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2037 {
2038   PetscErrorCode ierr;
2039   SNES           snes;
2040 
2041   PetscFunctionBegin;
2042   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2043   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2044   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2045   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2046   PetscFunctionReturn(0);
2047 }
2048 
2049 #undef __FUNCT__
2050 #define __FUNCT__ "TSGetOptionsPrefix"
2051 /*@C
2052    TSGetOptionsPrefix - Sets the prefix used for searching for all
2053    TS options in the database.
2054 
2055    Not Collective
2056 
2057    Input Parameter:
2058 .  ts - The TS context
2059 
2060    Output Parameter:
2061 .  prefix - A pointer to the prefix string used
2062 
2063    Notes: On the fortran side, the user should pass in a string 'prifix' of
2064    sufficient length to hold the prefix.
2065 
2066    Level: intermediate
2067 
2068 .keywords: TS, get, options, prefix, database
2069 
2070 .seealso: TSAppendOptionsPrefix()
2071 @*/
2072 PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2073 {
2074   PetscErrorCode ierr;
2075 
2076   PetscFunctionBegin;
2077   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2078   PetscValidPointer(prefix,2);
2079   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2080   PetscFunctionReturn(0);
2081 }
2082 
2083 #undef __FUNCT__
2084 #define __FUNCT__ "TSGetRHSJacobian"
2085 /*@C
2086    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2087 
2088    Not Collective, but parallel objects are returned if TS is parallel
2089 
2090    Input Parameter:
2091 .  ts  - The TS context obtained from TSCreate()
2092 
2093    Output Parameters:
2094 +  J   - The Jacobian J of F, where U_t = F(U,t)
2095 .  M   - The preconditioner matrix, usually the same as J
2096 .  func - Function to compute the Jacobian of the RHS
2097 -  ctx - User-defined context for Jacobian evaluation routine
2098 
2099    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2100 
2101    Level: intermediate
2102 
2103 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2104 
2105 .keywords: TS, timestep, get, matrix, Jacobian
2106 @*/
2107 PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2108 {
2109   PetscErrorCode ierr;
2110   SNES           snes;
2111 
2112   PetscFunctionBegin;
2113   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2114   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2115   if (func) *func = ts->ops->rhsjacobian;
2116   if (ctx) *ctx = ts->jacP;
2117   PetscFunctionReturn(0);
2118 }
2119 
2120 #undef __FUNCT__
2121 #define __FUNCT__ "TSGetIJacobian"
2122 /*@C
2123    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
2124 
2125    Not Collective, but parallel objects are returned if TS is parallel
2126 
2127    Input Parameter:
2128 .  ts  - The TS context obtained from TSCreate()
2129 
2130    Output Parameters:
2131 +  A   - The Jacobian of F(t,U,U_t)
2132 .  B   - The preconditioner matrix, often the same as A
2133 .  f   - The function to compute the matrices
2134 - ctx - User-defined context for Jacobian evaluation routine
2135 
2136    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2137 
2138    Level: advanced
2139 
2140 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2141 
2142 .keywords: TS, timestep, get, matrix, Jacobian
2143 @*/
2144 PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2145 {
2146   PetscErrorCode ierr;
2147   SNES           snes;
2148 
2149   PetscFunctionBegin;
2150   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2151   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2152   if (f) *f = ts->ops->ijacobian;
2153   if (ctx) *ctx = ts->jacP;
2154   PetscFunctionReturn(0);
2155 }
2156 
2157 #undef __FUNCT__
2158 #define __FUNCT__ "TSMonitorSolution"
2159 /*@C
2160    TSMonitorSolution - Monitors progress of the TS solvers by calling
2161    VecView() for the solution at each timestep
2162 
2163    Collective on TS
2164 
2165    Input Parameters:
2166 +  ts - the TS context
2167 .  step - current time-step
2168 .  ptime - current time
2169 -  dummy - either a viewer or PETSC_NULL
2170 
2171    Level: intermediate
2172 
2173 .keywords: TS,  vector, monitor, view
2174 
2175 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2176 @*/
2177 PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2178 {
2179   PetscErrorCode ierr;
2180   PetscViewer    viewer = (PetscViewer) dummy;
2181 
2182   PetscFunctionBegin;
2183   if (!dummy) {
2184     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2185   }
2186   ierr = VecView(x,viewer);CHKERRQ(ierr);
2187   PetscFunctionReturn(0);
2188 }
2189 
2190 
2191 #undef __FUNCT__
2192 #define __FUNCT__ "TSSetDM"
2193 /*@
2194    TSSetDM - Sets the DM that may be used by some preconditioners
2195 
2196    Logically Collective on TS and DM
2197 
2198    Input Parameters:
2199 +  ts - the preconditioner context
2200 -  dm - the dm
2201 
2202    Level: intermediate
2203 
2204 
2205 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2206 @*/
2207 PetscErrorCode  TSSetDM(TS ts,DM dm)
2208 {
2209   PetscErrorCode ierr;
2210   SNES           snes;
2211 
2212   PetscFunctionBegin;
2213   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2214   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
2215   ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
2216   ts->dm = dm;
2217   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2218   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
2219   PetscFunctionReturn(0);
2220 }
2221 
2222 #undef __FUNCT__
2223 #define __FUNCT__ "TSGetDM"
2224 /*@
2225    TSGetDM - Gets the DM that may be used by some preconditioners
2226 
2227    Not Collective
2228 
2229    Input Parameter:
2230 . ts - the preconditioner context
2231 
2232    Output Parameter:
2233 .  dm - the dm
2234 
2235    Level: intermediate
2236 
2237 
2238 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2239 @*/
2240 PetscErrorCode  TSGetDM(TS ts,DM *dm)
2241 {
2242   PetscFunctionBegin;
2243   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2244   *dm = ts->dm;
2245   PetscFunctionReturn(0);
2246 }
2247 
2248 #undef __FUNCT__
2249 #define __FUNCT__ "SNESTSFormFunction"
2250 /*@
2251    SNESTSFormFunction - Function to evaluate nonlinear residual
2252 
2253    Logically Collective on SNES
2254 
2255    Input Parameter:
2256 + snes - nonlinear solver
2257 . X - the current state at which to evaluate the residual
2258 - ctx - user context, must be a TS
2259 
2260    Output Parameter:
2261 . F - the nonlinear residual
2262 
2263    Notes:
2264    This function is not normally called by users and is automatically registered with the SNES used by TS.
2265    It is most frequently passed to MatFDColoringSetFunction().
2266 
2267    Level: advanced
2268 
2269 .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2270 @*/
2271 PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2272 {
2273   TS ts = (TS)ctx;
2274   PetscErrorCode ierr;
2275 
2276   PetscFunctionBegin;
2277   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2278   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2279   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
2280   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
2281   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
2282   PetscFunctionReturn(0);
2283 }
2284 
2285 #undef __FUNCT__
2286 #define __FUNCT__ "SNESTSFormJacobian"
2287 /*@
2288    SNESTSFormJacobian - Function to evaluate the Jacobian
2289 
2290    Collective on SNES
2291 
2292    Input Parameter:
2293 + snes - nonlinear solver
2294 . X - the current state at which to evaluate the residual
2295 - ctx - user context, must be a TS
2296 
2297    Output Parameter:
2298 + A - the Jacobian
2299 . B - the preconditioning matrix (may be the same as A)
2300 - flag - indicates any structure change in the matrix
2301 
2302    Notes:
2303    This function is not normally called by users and is automatically registered with the SNES used by TS.
2304 
2305    Level: developer
2306 
2307 .seealso: SNESSetJacobian()
2308 @*/
2309 PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2310 {
2311   TS ts = (TS)ctx;
2312   PetscErrorCode ierr;
2313 
2314   PetscFunctionBegin;
2315   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2316   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2317   PetscValidPointer(A,3);
2318   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
2319   PetscValidPointer(B,4);
2320   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
2321   PetscValidPointer(flag,5);
2322   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
2323   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
2324   PetscFunctionReturn(0);
2325 }
2326 
2327 #if defined(PETSC_HAVE_MATLAB_ENGINE)
2328 #include <mex.h>
2329 
2330 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
2331 
2332 #undef __FUNCT__
2333 #define __FUNCT__ "TSComputeFunction_Matlab"
2334 /*
2335    TSComputeFunction_Matlab - Calls the function that has been set with
2336                          TSSetFunctionMatlab().
2337 
2338    Collective on TS
2339 
2340    Input Parameters:
2341 +  snes - the TS context
2342 -  x - input vector
2343 
2344    Output Parameter:
2345 .  y - function vector, as set by TSSetFunction()
2346 
2347    Notes:
2348    TSComputeFunction() is typically used within nonlinear solvers
2349    implementations, so most users would not generally call this routine
2350    themselves.
2351 
2352    Level: developer
2353 
2354 .keywords: TS, nonlinear, compute, function
2355 
2356 .seealso: TSSetFunction(), TSGetFunction()
2357 */
2358 PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
2359 {
2360   PetscErrorCode   ierr;
2361   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2362   int              nlhs = 1,nrhs = 7;
2363   mxArray	   *plhs[1],*prhs[7];
2364   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
2365 
2366   PetscFunctionBegin;
2367   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2368   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2369   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
2370   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
2371   PetscCheckSameComm(snes,1,x,3);
2372   PetscCheckSameComm(snes,1,y,5);
2373 
2374   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2375   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2376   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
2377   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
2378   prhs[0] =  mxCreateDoubleScalar((double)ls);
2379   prhs[1] =  mxCreateDoubleScalar(time);
2380   prhs[2] =  mxCreateDoubleScalar((double)lx);
2381   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2382   prhs[4] =  mxCreateDoubleScalar((double)ly);
2383   prhs[5] =  mxCreateString(sctx->funcname);
2384   prhs[6] =  sctx->ctx;
2385   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
2386   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2387   mxDestroyArray(prhs[0]);
2388   mxDestroyArray(prhs[1]);
2389   mxDestroyArray(prhs[2]);
2390   mxDestroyArray(prhs[3]);
2391   mxDestroyArray(prhs[4]);
2392   mxDestroyArray(prhs[5]);
2393   mxDestroyArray(plhs[0]);
2394   PetscFunctionReturn(0);
2395 }
2396 
2397 
2398 #undef __FUNCT__
2399 #define __FUNCT__ "TSSetFunctionMatlab"
2400 /*
2401    TSSetFunctionMatlab - Sets the function evaluation routine and function
2402    vector for use by the TS routines in solving ODEs
2403    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
2404 
2405    Logically Collective on TS
2406 
2407    Input Parameters:
2408 +  ts - the TS context
2409 -  func - function evaluation routine
2410 
2411    Calling sequence of func:
2412 $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
2413 
2414    Level: beginner
2415 
2416 .keywords: TS, nonlinear, set, function
2417 
2418 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2419 */
2420 PetscErrorCode  TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx)
2421 {
2422   PetscErrorCode  ierr;
2423   TSMatlabContext *sctx;
2424 
2425   PetscFunctionBegin;
2426   /* currently sctx is memory bleed */
2427   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2428   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2429   /*
2430      This should work, but it doesn't
2431   sctx->ctx = ctx;
2432   mexMakeArrayPersistent(sctx->ctx);
2433   */
2434   sctx->ctx = mxDuplicateArray(ctx);
2435   ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
2436   PetscFunctionReturn(0);
2437 }
2438 
2439 #undef __FUNCT__
2440 #define __FUNCT__ "TSComputeJacobian_Matlab"
2441 /*
2442    TSComputeJacobian_Matlab - Calls the function that has been set with
2443                          TSSetJacobianMatlab().
2444 
2445    Collective on TS
2446 
2447    Input Parameters:
2448 +  snes - the TS context
2449 .  x - input vector
2450 .  A, B - the matrices
2451 -  ctx - user context
2452 
2453    Output Parameter:
2454 .  flag - structure of the matrix
2455 
2456    Level: developer
2457 
2458 .keywords: TS, nonlinear, compute, function
2459 
2460 .seealso: TSSetFunction(), TSGetFunction()
2461 @*/
2462 PetscErrorCode  TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
2463 {
2464   PetscErrorCode  ierr;
2465   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2466   int             nlhs = 2,nrhs = 9;
2467   mxArray	  *plhs[2],*prhs[9];
2468   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
2469 
2470   PetscFunctionBegin;
2471   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2472   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2473 
2474   /* call Matlab function in ctx with arguments x and y */
2475 
2476   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2477   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2478   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
2479   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
2480   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
2481   prhs[0] =  mxCreateDoubleScalar((double)ls);
2482   prhs[1] =  mxCreateDoubleScalar((double)time);
2483   prhs[2] =  mxCreateDoubleScalar((double)lx);
2484   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2485   prhs[4] =  mxCreateDoubleScalar((double)shift);
2486   prhs[5] =  mxCreateDoubleScalar((double)lA);
2487   prhs[6] =  mxCreateDoubleScalar((double)lB);
2488   prhs[7] =  mxCreateString(sctx->funcname);
2489   prhs[8] =  sctx->ctx;
2490   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
2491   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2492   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
2493   mxDestroyArray(prhs[0]);
2494   mxDestroyArray(prhs[1]);
2495   mxDestroyArray(prhs[2]);
2496   mxDestroyArray(prhs[3]);
2497   mxDestroyArray(prhs[4]);
2498   mxDestroyArray(prhs[5]);
2499   mxDestroyArray(prhs[6]);
2500   mxDestroyArray(prhs[7]);
2501   mxDestroyArray(plhs[0]);
2502   mxDestroyArray(plhs[1]);
2503   PetscFunctionReturn(0);
2504 }
2505 
2506 
2507 #undef __FUNCT__
2508 #define __FUNCT__ "TSSetJacobianMatlab"
2509 /*
2510    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
2511    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
2512 
2513    Logically Collective on TS
2514 
2515    Input Parameters:
2516 +  snes - the TS context
2517 .  A,B - Jacobian matrices
2518 .  func - function evaluation routine
2519 -  ctx - user context
2520 
2521    Calling sequence of func:
2522 $    flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
2523 
2524 
2525    Level: developer
2526 
2527 .keywords: TS, nonlinear, set, function
2528 
2529 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2530 */
2531 PetscErrorCode  TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx)
2532 {
2533   PetscErrorCode    ierr;
2534   TSMatlabContext *sctx;
2535 
2536   PetscFunctionBegin;
2537   /* currently sctx is memory bleed */
2538   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2539   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2540   /*
2541      This should work, but it doesn't
2542   sctx->ctx = ctx;
2543   mexMakeArrayPersistent(sctx->ctx);
2544   */
2545   sctx->ctx = mxDuplicateArray(ctx);
2546   ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
2547   PetscFunctionReturn(0);
2548 }
2549 
2550 #undef __FUNCT__
2551 #define __FUNCT__ "TSMonitor_Matlab"
2552 /*
2553    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
2554 
2555    Collective on TS
2556 
2557 .seealso: TSSetFunction(), TSGetFunction()
2558 @*/
2559 PetscErrorCode  TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx)
2560 {
2561   PetscErrorCode  ierr;
2562   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2563   int             nlhs = 1,nrhs = 6;
2564   mxArray	  *plhs[1],*prhs[6];
2565   long long int   lx = 0,ls = 0;
2566 
2567   PetscFunctionBegin;
2568   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2569   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
2570 
2571   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2572   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2573   prhs[0] =  mxCreateDoubleScalar((double)ls);
2574   prhs[1] =  mxCreateDoubleScalar((double)it);
2575   prhs[2] =  mxCreateDoubleScalar((double)time);
2576   prhs[3] =  mxCreateDoubleScalar((double)lx);
2577   prhs[4] =  mxCreateString(sctx->funcname);
2578   prhs[5] =  sctx->ctx;
2579   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
2580   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2581   mxDestroyArray(prhs[0]);
2582   mxDestroyArray(prhs[1]);
2583   mxDestroyArray(prhs[2]);
2584   mxDestroyArray(prhs[3]);
2585   mxDestroyArray(prhs[4]);
2586   mxDestroyArray(plhs[0]);
2587   PetscFunctionReturn(0);
2588 }
2589 
2590 
2591 #undef __FUNCT__
2592 #define __FUNCT__ "TSMonitorSetMatlab"
2593 /*
2594    TSMonitorSetMatlab - Sets the monitor function from Matlab
2595 
2596    Level: developer
2597 
2598 .keywords: TS, nonlinear, set, function
2599 
2600 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2601 */
2602 PetscErrorCode  TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx)
2603 {
2604   PetscErrorCode    ierr;
2605   TSMatlabContext *sctx;
2606 
2607   PetscFunctionBegin;
2608   /* currently sctx is memory bleed */
2609   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2610   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2611   /*
2612      This should work, but it doesn't
2613   sctx->ctx = ctx;
2614   mexMakeArrayPersistent(sctx->ctx);
2615   */
2616   sctx->ctx = mxDuplicateArray(ctx);
2617   ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
2618   PetscFunctionReturn(0);
2619 }
2620 
2621 #endif
2622