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