xref: /petsc/src/ts/interface/ts.c (revision cd6526761e44e2ac7b4afaaf511785e78693abcc)
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__ "TSSetRetainStages"
1650 /*@
1651    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
1652 
1653    Logically Collective on TS
1654 
1655    Input Argument:
1656 .  ts - time stepping context
1657 
1658    Output Argument:
1659 .  flg - PETSC_TRUE or PETSC_FALSE
1660 
1661    Level: intermediate
1662 
1663 .keywords: TS, set
1664 
1665 .seealso: TSInterpolate(), TSSetPostStep()
1666 @*/
1667 PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
1668 {
1669 
1670   PetscFunctionBegin;
1671   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1672   ts->retain_stages = flg;
1673   PetscFunctionReturn(0);
1674 }
1675 
1676 #undef __FUNCT__
1677 #define __FUNCT__ "TSInterpolate"
1678 /*@
1679    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
1680 
1681    Collective on TS
1682 
1683    Input Argument:
1684 +  ts - time stepping context
1685 -  t - time to interpolate to
1686 
1687    Output Argument:
1688 .  X - state at given time
1689 
1690    Notes:
1691    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
1692 
1693    Level: intermediate
1694 
1695    Developer Notes:
1696    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
1697 
1698 .keywords: TS, set
1699 
1700 .seealso: TSSetRetainStages(), TSSetPostStep()
1701 @*/
1702 PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
1703 {
1704   PetscErrorCode ierr;
1705 
1706   PetscFunctionBegin;
1707   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1708   if (t < ts->ptime - ts->time_step || ts->ptime < t) SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested time %G not in last time steps [%G,%G]",t,ts->ptime-ts->time_step,ts->ptime);
1709   if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
1710   ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr);
1711   PetscFunctionReturn(0);
1712 }
1713 
1714 #undef __FUNCT__
1715 #define __FUNCT__ "TSStep"
1716 /*@
1717    TSStep - Steps the requested number of timesteps.
1718 
1719    Collective on TS
1720 
1721    Input Parameter:
1722 .  ts - the TS context obtained from TSCreate()
1723 
1724    Output Parameters:
1725 +  steps - number of iterations until termination
1726 -  ptime - time until termination
1727 
1728    Level: beginner
1729 
1730 .keywords: TS, timestep, solve
1731 
1732 .seealso: TSCreate(), TSSetUp(), TSDestroy()
1733 @*/
1734 PetscErrorCode  TSStep(TS ts)
1735 {
1736   PetscErrorCode ierr;
1737 
1738   PetscFunctionBegin;
1739   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1740 
1741   ierr = TSSetUp(ts);CHKERRQ(ierr);
1742 
1743   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1744   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
1745   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1746   PetscFunctionReturn(0);
1747 }
1748 
1749 #undef __FUNCT__
1750 #define __FUNCT__ "TSSolve"
1751 /*@
1752    TSSolve - Steps the requested number of timesteps.
1753 
1754    Collective on TS
1755 
1756    Input Parameter:
1757 +  ts - the TS context obtained from TSCreate()
1758 -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
1759 
1760    Level: beginner
1761 
1762 .keywords: TS, timestep, solve
1763 
1764 .seealso: TSCreate(), TSSetSolution(), TSStep()
1765 @*/
1766 PetscErrorCode  TSSolve(TS ts, Vec x)
1767 {
1768   PetscInt       i;
1769   PetscBool      flg;
1770   char           filename[PETSC_MAX_PATH_LEN];
1771   PetscViewer    viewer;
1772   PetscErrorCode ierr;
1773 
1774   PetscFunctionBegin;
1775   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1776   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1777   ierr = TSSetSolution(ts,x); CHKERRQ(ierr);
1778   ierr = TSSetUp(ts); CHKERRQ(ierr);
1779   /* reset time step and iteration counters */
1780   ts->steps = 0;
1781   ts->linear_its = 0;
1782   ts->nonlinear_its = 0;
1783   ts->reason = TS_CONVERGED_ITERATING;
1784   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1785 
1786   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
1787     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
1788   } else {
1789     i = 0;
1790     if (i >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
1791     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
1792     /* steps the requested number of timesteps. */
1793     while (!ts->reason) {
1794       ierr = TSPreStep(ts);CHKERRQ(ierr);
1795       ierr = TSStep(ts);CHKERRQ(ierr);
1796       if (ts->reason < 0) {
1797         if (ts->errorifstepfailed) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed");
1798       } else if (++i >= ts->max_steps) {
1799         ts->reason = TS_CONVERGED_ITS;
1800       } else if (ts->ptime >= ts->max_time) {
1801         ts->reason = TS_CONVERGED_TIME;
1802       }
1803       ierr = TSPostStep(ts);CHKERRQ(ierr);
1804       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1805     }
1806   }
1807   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1808   if (flg && !PetscPreLoadingOn) {
1809     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr);
1810     ierr = TSView(ts,viewer);CHKERRQ(ierr);
1811     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1812   }
1813   PetscFunctionReturn(0);
1814 }
1815 
1816 #undef __FUNCT__
1817 #define __FUNCT__ "TSMonitor"
1818 /*
1819      Runs the user provided monitor routines, if they exists.
1820 */
1821 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1822 {
1823   PetscErrorCode ierr;
1824   PetscInt       i,n = ts->numbermonitors;
1825 
1826   PetscFunctionBegin;
1827   for (i=0; i<n; i++) {
1828     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1829   }
1830   PetscFunctionReturn(0);
1831 }
1832 
1833 /* ------------------------------------------------------------------------*/
1834 
1835 #undef __FUNCT__
1836 #define __FUNCT__ "TSMonitorLGCreate"
1837 /*@C
1838    TSMonitorLGCreate - Creates a line graph context for use with
1839    TS to monitor convergence of preconditioned residual norms.
1840 
1841    Collective on TS
1842 
1843    Input Parameters:
1844 +  host - the X display to open, or null for the local machine
1845 .  label - the title to put in the title bar
1846 .  x, y - the screen coordinates of the upper left coordinate of the window
1847 -  m, n - the screen width and height in pixels
1848 
1849    Output Parameter:
1850 .  draw - the drawing context
1851 
1852    Options Database Key:
1853 .  -ts_monitor_draw - automatically sets line graph monitor
1854 
1855    Notes:
1856    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1857 
1858    Level: intermediate
1859 
1860 .keywords: TS, monitor, line graph, residual, seealso
1861 
1862 .seealso: TSMonitorLGDestroy(), TSMonitorSet()
1863 
1864 @*/
1865 PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1866 {
1867   PetscDraw      win;
1868   PetscErrorCode ierr;
1869 
1870   PetscFunctionBegin;
1871   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1872   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1873   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1874   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1875 
1876   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1877   PetscFunctionReturn(0);
1878 }
1879 
1880 #undef __FUNCT__
1881 #define __FUNCT__ "TSMonitorLG"
1882 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1883 {
1884   PetscDrawLG    lg = (PetscDrawLG) monctx;
1885   PetscReal      x,y = ptime;
1886   PetscErrorCode ierr;
1887 
1888   PetscFunctionBegin;
1889   if (!monctx) {
1890     MPI_Comm    comm;
1891     PetscViewer viewer;
1892 
1893     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1894     viewer = PETSC_VIEWER_DRAW_(comm);
1895     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
1896   }
1897 
1898   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1899   x = (PetscReal)n;
1900   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1901   if (n < 20 || (n % 5)) {
1902     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1903   }
1904   PetscFunctionReturn(0);
1905 }
1906 
1907 #undef __FUNCT__
1908 #define __FUNCT__ "TSMonitorLGDestroy"
1909 /*@C
1910    TSMonitorLGDestroy - Destroys a line graph context that was created
1911    with TSMonitorLGCreate().
1912 
1913    Collective on PetscDrawLG
1914 
1915    Input Parameter:
1916 .  draw - the drawing context
1917 
1918    Level: intermediate
1919 
1920 .keywords: TS, monitor, line graph, destroy
1921 
1922 .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1923 @*/
1924 PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
1925 {
1926   PetscDraw      draw;
1927   PetscErrorCode ierr;
1928 
1929   PetscFunctionBegin;
1930   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
1931   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
1932   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1933   PetscFunctionReturn(0);
1934 }
1935 
1936 #undef __FUNCT__
1937 #define __FUNCT__ "TSGetTime"
1938 /*@
1939    TSGetTime - Gets the current time.
1940 
1941    Not Collective
1942 
1943    Input Parameter:
1944 .  ts - the TS context obtained from TSCreate()
1945 
1946    Output Parameter:
1947 .  t  - the current time
1948 
1949    Level: beginner
1950 
1951 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1952 
1953 .keywords: TS, get, time
1954 @*/
1955 PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
1956 {
1957   PetscFunctionBegin;
1958   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1959   PetscValidDoublePointer(t,2);
1960   *t = ts->ptime;
1961   PetscFunctionReturn(0);
1962 }
1963 
1964 #undef __FUNCT__
1965 #define __FUNCT__ "TSSetTime"
1966 /*@
1967    TSSetTime - Allows one to reset the time.
1968 
1969    Logically Collective on TS
1970 
1971    Input Parameters:
1972 +  ts - the TS context obtained from TSCreate()
1973 -  time - the time
1974 
1975    Level: intermediate
1976 
1977 .seealso: TSGetTime(), TSSetDuration()
1978 
1979 .keywords: TS, set, time
1980 @*/
1981 PetscErrorCode  TSSetTime(TS ts, PetscReal t)
1982 {
1983   PetscFunctionBegin;
1984   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1985   PetscValidLogicalCollectiveReal(ts,t,2);
1986   ts->ptime = t;
1987   PetscFunctionReturn(0);
1988 }
1989 
1990 #undef __FUNCT__
1991 #define __FUNCT__ "TSSetOptionsPrefix"
1992 /*@C
1993    TSSetOptionsPrefix - Sets the prefix used for searching for all
1994    TS options in the database.
1995 
1996    Logically Collective on TS
1997 
1998    Input Parameter:
1999 +  ts     - The TS context
2000 -  prefix - The prefix to prepend to all option names
2001 
2002    Notes:
2003    A hyphen (-) must NOT be given at the beginning of the prefix name.
2004    The first character of all runtime options is AUTOMATICALLY the
2005    hyphen.
2006 
2007    Level: advanced
2008 
2009 .keywords: TS, set, options, prefix, database
2010 
2011 .seealso: TSSetFromOptions()
2012 
2013 @*/
2014 PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2015 {
2016   PetscErrorCode ierr;
2017   SNES           snes;
2018 
2019   PetscFunctionBegin;
2020   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2021   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2022   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2023   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2024   PetscFunctionReturn(0);
2025 }
2026 
2027 
2028 #undef __FUNCT__
2029 #define __FUNCT__ "TSAppendOptionsPrefix"
2030 /*@C
2031    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2032    TS options in the database.
2033 
2034    Logically Collective on TS
2035 
2036    Input Parameter:
2037 +  ts     - The TS context
2038 -  prefix - The prefix to prepend to all option names
2039 
2040    Notes:
2041    A hyphen (-) must NOT be given at the beginning of the prefix name.
2042    The first character of all runtime options is AUTOMATICALLY the
2043    hyphen.
2044 
2045    Level: advanced
2046 
2047 .keywords: TS, append, options, prefix, database
2048 
2049 .seealso: TSGetOptionsPrefix()
2050 
2051 @*/
2052 PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2053 {
2054   PetscErrorCode ierr;
2055   SNES           snes;
2056 
2057   PetscFunctionBegin;
2058   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2059   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2060   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2061   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2062   PetscFunctionReturn(0);
2063 }
2064 
2065 #undef __FUNCT__
2066 #define __FUNCT__ "TSGetOptionsPrefix"
2067 /*@C
2068    TSGetOptionsPrefix - Sets the prefix used for searching for all
2069    TS options in the database.
2070 
2071    Not Collective
2072 
2073    Input Parameter:
2074 .  ts - The TS context
2075 
2076    Output Parameter:
2077 .  prefix - A pointer to the prefix string used
2078 
2079    Notes: On the fortran side, the user should pass in a string 'prifix' of
2080    sufficient length to hold the prefix.
2081 
2082    Level: intermediate
2083 
2084 .keywords: TS, get, options, prefix, database
2085 
2086 .seealso: TSAppendOptionsPrefix()
2087 @*/
2088 PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2089 {
2090   PetscErrorCode ierr;
2091 
2092   PetscFunctionBegin;
2093   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2094   PetscValidPointer(prefix,2);
2095   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2096   PetscFunctionReturn(0);
2097 }
2098 
2099 #undef __FUNCT__
2100 #define __FUNCT__ "TSGetRHSJacobian"
2101 /*@C
2102    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2103 
2104    Not Collective, but parallel objects are returned if TS is parallel
2105 
2106    Input Parameter:
2107 .  ts  - The TS context obtained from TSCreate()
2108 
2109    Output Parameters:
2110 +  J   - The Jacobian J of F, where U_t = F(U,t)
2111 .  M   - The preconditioner matrix, usually the same as J
2112 .  func - Function to compute the Jacobian of the RHS
2113 -  ctx - User-defined context for Jacobian evaluation routine
2114 
2115    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2116 
2117    Level: intermediate
2118 
2119 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2120 
2121 .keywords: TS, timestep, get, matrix, Jacobian
2122 @*/
2123 PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2124 {
2125   PetscErrorCode ierr;
2126   SNES           snes;
2127 
2128   PetscFunctionBegin;
2129   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2130   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2131   if (func) *func = ts->userops->rhsjacobian;
2132   if (ctx) *ctx = ts->jacP;
2133   PetscFunctionReturn(0);
2134 }
2135 
2136 #undef __FUNCT__
2137 #define __FUNCT__ "TSGetIJacobian"
2138 /*@C
2139    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
2140 
2141    Not Collective, but parallel objects are returned if TS is parallel
2142 
2143    Input Parameter:
2144 .  ts  - The TS context obtained from TSCreate()
2145 
2146    Output Parameters:
2147 +  A   - The Jacobian of F(t,U,U_t)
2148 .  B   - The preconditioner matrix, often the same as A
2149 .  f   - The function to compute the matrices
2150 - ctx - User-defined context for Jacobian evaluation routine
2151 
2152    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2153 
2154    Level: advanced
2155 
2156 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2157 
2158 .keywords: TS, timestep, get, matrix, Jacobian
2159 @*/
2160 PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2161 {
2162   PetscErrorCode ierr;
2163   SNES           snes;
2164 
2165   PetscFunctionBegin;
2166   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2167   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2168   if (f) *f = ts->userops->ijacobian;
2169   if (ctx) *ctx = ts->jacP;
2170   PetscFunctionReturn(0);
2171 }
2172 
2173 #undef __FUNCT__
2174 #define __FUNCT__ "TSMonitorSolution"
2175 /*@C
2176    TSMonitorSolution - Monitors progress of the TS solvers by calling
2177    VecView() for the solution at each timestep
2178 
2179    Collective on TS
2180 
2181    Input Parameters:
2182 +  ts - the TS context
2183 .  step - current time-step
2184 .  ptime - current time
2185 -  dummy - either a viewer or PETSC_NULL
2186 
2187    Level: intermediate
2188 
2189 .keywords: TS,  vector, monitor, view
2190 
2191 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2192 @*/
2193 PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2194 {
2195   PetscErrorCode ierr;
2196   PetscViewer    viewer = (PetscViewer) dummy;
2197 
2198   PetscFunctionBegin;
2199   if (!dummy) {
2200     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2201   }
2202   ierr = VecView(x,viewer);CHKERRQ(ierr);
2203   PetscFunctionReturn(0);
2204 }
2205 
2206 
2207 #undef __FUNCT__
2208 #define __FUNCT__ "TSSetDM"
2209 /*@
2210    TSSetDM - Sets the DM that may be used by some preconditioners
2211 
2212    Logically Collective on TS and DM
2213 
2214    Input Parameters:
2215 +  ts - the preconditioner context
2216 -  dm - the dm
2217 
2218    Level: intermediate
2219 
2220 
2221 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2222 @*/
2223 PetscErrorCode  TSSetDM(TS ts,DM dm)
2224 {
2225   PetscErrorCode ierr;
2226   SNES           snes;
2227 
2228   PetscFunctionBegin;
2229   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2230   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
2231   ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
2232   ts->dm = dm;
2233   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2234   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
2235   PetscFunctionReturn(0);
2236 }
2237 
2238 #undef __FUNCT__
2239 #define __FUNCT__ "TSGetDM"
2240 /*@
2241    TSGetDM - Gets the DM that may be used by some preconditioners
2242 
2243    Not Collective
2244 
2245    Input Parameter:
2246 . ts - the preconditioner context
2247 
2248    Output Parameter:
2249 .  dm - the dm
2250 
2251    Level: intermediate
2252 
2253 
2254 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2255 @*/
2256 PetscErrorCode  TSGetDM(TS ts,DM *dm)
2257 {
2258   PetscFunctionBegin;
2259   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2260   *dm = ts->dm;
2261   PetscFunctionReturn(0);
2262 }
2263 
2264 #undef __FUNCT__
2265 #define __FUNCT__ "SNESTSFormFunction"
2266 /*@
2267    SNESTSFormFunction - Function to evaluate nonlinear residual
2268 
2269    Logically Collective on SNES
2270 
2271    Input Parameter:
2272 + snes - nonlinear solver
2273 . X - the current state at which to evaluate the residual
2274 - ctx - user context, must be a TS
2275 
2276    Output Parameter:
2277 . F - the nonlinear residual
2278 
2279    Notes:
2280    This function is not normally called by users and is automatically registered with the SNES used by TS.
2281    It is most frequently passed to MatFDColoringSetFunction().
2282 
2283    Level: advanced
2284 
2285 .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2286 @*/
2287 PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2288 {
2289   TS ts = (TS)ctx;
2290   PetscErrorCode ierr;
2291 
2292   PetscFunctionBegin;
2293   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2294   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2295   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
2296   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
2297   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
2298   PetscFunctionReturn(0);
2299 }
2300 
2301 #undef __FUNCT__
2302 #define __FUNCT__ "SNESTSFormJacobian"
2303 /*@
2304    SNESTSFormJacobian - Function to evaluate the Jacobian
2305 
2306    Collective on SNES
2307 
2308    Input Parameter:
2309 + snes - nonlinear solver
2310 . X - the current state at which to evaluate the residual
2311 - ctx - user context, must be a TS
2312 
2313    Output Parameter:
2314 + A - the Jacobian
2315 . B - the preconditioning matrix (may be the same as A)
2316 - flag - indicates any structure change in the matrix
2317 
2318    Notes:
2319    This function is not normally called by users and is automatically registered with the SNES used by TS.
2320 
2321    Level: developer
2322 
2323 .seealso: SNESSetJacobian()
2324 @*/
2325 PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2326 {
2327   TS ts = (TS)ctx;
2328   PetscErrorCode ierr;
2329 
2330   PetscFunctionBegin;
2331   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2332   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2333   PetscValidPointer(A,3);
2334   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
2335   PetscValidPointer(B,4);
2336   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
2337   PetscValidPointer(flag,5);
2338   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
2339   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
2340   PetscFunctionReturn(0);
2341 }
2342 
2343 #undef __FUNCT__
2344 #define __FUNCT__ "TSComputeRHSFunctionLinear"
2345 /*@C
2346    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only
2347 
2348    Collective on TS
2349 
2350    Input Arguments:
2351 +  ts - time stepping context
2352 .  t - time at which to evaluate
2353 .  X - state at which to evaluate
2354 -  ctx - context
2355 
2356    Output Arguments:
2357 .  F - right hand side
2358 
2359    Level: intermediate
2360 
2361    Notes:
2362    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
2363    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
2364 
2365 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
2366 @*/
2367 PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
2368 {
2369   PetscErrorCode ierr;
2370   Mat Arhs,Brhs;
2371   MatStructure flg2;
2372 
2373   PetscFunctionBegin;
2374   ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
2375   ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
2376   ierr = MatMult(Arhs,X,F);CHKERRQ(ierr);
2377   PetscFunctionReturn(0);
2378 }
2379 
2380 #undef __FUNCT__
2381 #define __FUNCT__ "TSComputeRHSJacobianConstant"
2382 /*@C
2383    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
2384 
2385    Collective on TS
2386 
2387    Input Arguments:
2388 +  ts - time stepping context
2389 .  t - time at which to evaluate
2390 .  X - state at which to evaluate
2391 -  ctx - context
2392 
2393    Output Arguments:
2394 +  A - pointer to operator
2395 .  B - pointer to preconditioning matrix
2396 -  flg - matrix structure flag
2397 
2398    Level: intermediate
2399 
2400    Notes:
2401    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
2402 
2403 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
2404 @*/
2405 PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2406 {
2407 
2408   PetscFunctionBegin;
2409   *flg = SAME_PRECONDITIONER;
2410   PetscFunctionReturn(0);
2411 }
2412 
2413 #undef __FUNCT__
2414 #define __FUNCT__ "TSComputeIFunctionLinear"
2415 /*@C
2416    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
2417 
2418    Collective on TS
2419 
2420    Input Arguments:
2421 +  ts - time stepping context
2422 .  t - time at which to evaluate
2423 .  X - state at which to evaluate
2424 .  Xdot - time derivative of state vector
2425 -  ctx - context
2426 
2427    Output Arguments:
2428 .  F - left hand side
2429 
2430    Level: intermediate
2431 
2432    Notes:
2433    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
2434    user is required to write their own TSComputeIFunction.
2435    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
2436    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
2437 
2438 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
2439 @*/
2440 PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
2441 {
2442   PetscErrorCode ierr;
2443   Mat A,B;
2444   MatStructure flg2;
2445 
2446   PetscFunctionBegin;
2447   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2448   ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr);
2449   ierr = MatMult(A,Xdot,F);CHKERRQ(ierr);
2450   PetscFunctionReturn(0);
2451 }
2452 
2453 #undef __FUNCT__
2454 #define __FUNCT__ "TSComputeIJacobianConstant"
2455 /*@C
2456    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
2457 
2458    Collective on TS
2459 
2460    Input Arguments:
2461 +  ts - time stepping context
2462 .  t - time at which to evaluate
2463 .  X - state at which to evaluate
2464 .  Xdot - time derivative of state vector
2465 .  shift - shift to apply
2466 -  ctx - context
2467 
2468    Output Arguments:
2469 +  A - pointer to operator
2470 .  B - pointer to preconditioning matrix
2471 -  flg - matrix structure flag
2472 
2473    Level: intermediate
2474 
2475    Notes:
2476    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
2477 
2478 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
2479 @*/
2480 PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2481 {
2482 
2483   PetscFunctionBegin;
2484   *flg = SAME_PRECONDITIONER;
2485   PetscFunctionReturn(0);
2486 }
2487 
2488 
2489 #undef __FUNCT__
2490 #define __FUNCT__ "TSGetConvergedReason"
2491 /*@
2492    TSGetConvergedReason - Gets the reason the TS iteration was stopped.
2493 
2494    Not Collective
2495 
2496    Input Parameter:
2497 .  ts - the TS context
2498 
2499    Output Parameter:
2500 .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
2501             manual pages for the individual convergence tests for complete lists
2502 
2503    Level: intermediate
2504 
2505    Notes:
2506    Can only be called after the call to TSSolve() is complete.
2507 
2508 .keywords: TS, nonlinear, set, convergence, test
2509 
2510 .seealso: TSSetConvergenceTest(), TSConvergedReason
2511 @*/
2512 PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
2513 {
2514   PetscFunctionBegin;
2515   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2516   PetscValidPointer(reason,2);
2517   *reason = ts->reason;
2518   PetscFunctionReturn(0);
2519 }
2520 
2521 #if defined(PETSC_HAVE_MATLAB_ENGINE)
2522 #include <mex.h>
2523 
2524 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
2525 
2526 #undef __FUNCT__
2527 #define __FUNCT__ "TSComputeFunction_Matlab"
2528 /*
2529    TSComputeFunction_Matlab - Calls the function that has been set with
2530                          TSSetFunctionMatlab().
2531 
2532    Collective on TS
2533 
2534    Input Parameters:
2535 +  snes - the TS context
2536 -  x - input vector
2537 
2538    Output Parameter:
2539 .  y - function vector, as set by TSSetFunction()
2540 
2541    Notes:
2542    TSComputeFunction() is typically used within nonlinear solvers
2543    implementations, so most users would not generally call this routine
2544    themselves.
2545 
2546    Level: developer
2547 
2548 .keywords: TS, nonlinear, compute, function
2549 
2550 .seealso: TSSetFunction(), TSGetFunction()
2551 */
2552 PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
2553 {
2554   PetscErrorCode   ierr;
2555   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2556   int              nlhs = 1,nrhs = 7;
2557   mxArray          *plhs[1],*prhs[7];
2558   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
2559 
2560   PetscFunctionBegin;
2561   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2562   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2563   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
2564   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
2565   PetscCheckSameComm(snes,1,x,3);
2566   PetscCheckSameComm(snes,1,y,5);
2567 
2568   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2569   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2570   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
2571   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
2572   prhs[0] =  mxCreateDoubleScalar((double)ls);
2573   prhs[1] =  mxCreateDoubleScalar(time);
2574   prhs[2] =  mxCreateDoubleScalar((double)lx);
2575   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2576   prhs[4] =  mxCreateDoubleScalar((double)ly);
2577   prhs[5] =  mxCreateString(sctx->funcname);
2578   prhs[6] =  sctx->ctx;
2579   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
2580   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2581   mxDestroyArray(prhs[0]);
2582   mxDestroyArray(prhs[1]);
2583   mxDestroyArray(prhs[2]);
2584   mxDestroyArray(prhs[3]);
2585   mxDestroyArray(prhs[4]);
2586   mxDestroyArray(prhs[5]);
2587   mxDestroyArray(plhs[0]);
2588   PetscFunctionReturn(0);
2589 }
2590 
2591 
2592 #undef __FUNCT__
2593 #define __FUNCT__ "TSSetFunctionMatlab"
2594 /*
2595    TSSetFunctionMatlab - Sets the function evaluation routine and function
2596    vector for use by the TS routines in solving ODEs
2597    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
2598 
2599    Logically Collective on TS
2600 
2601    Input Parameters:
2602 +  ts - the TS context
2603 -  func - function evaluation routine
2604 
2605    Calling sequence of func:
2606 $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
2607 
2608    Level: beginner
2609 
2610 .keywords: TS, nonlinear, set, function
2611 
2612 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2613 */
2614 PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
2615 {
2616   PetscErrorCode  ierr;
2617   TSMatlabContext *sctx;
2618 
2619   PetscFunctionBegin;
2620   /* currently sctx is memory bleed */
2621   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2622   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2623   /*
2624      This should work, but it doesn't
2625   sctx->ctx = ctx;
2626   mexMakeArrayPersistent(sctx->ctx);
2627   */
2628   sctx->ctx = mxDuplicateArray(ctx);
2629   ierr = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
2630   PetscFunctionReturn(0);
2631 }
2632 
2633 #undef __FUNCT__
2634 #define __FUNCT__ "TSComputeJacobian_Matlab"
2635 /*
2636    TSComputeJacobian_Matlab - Calls the function that has been set with
2637                          TSSetJacobianMatlab().
2638 
2639    Collective on TS
2640 
2641    Input Parameters:
2642 +  ts - the TS context
2643 .  x - input vector
2644 .  A, B - the matrices
2645 -  ctx - user context
2646 
2647    Output Parameter:
2648 .  flag - structure of the matrix
2649 
2650    Level: developer
2651 
2652 .keywords: TS, nonlinear, compute, function
2653 
2654 .seealso: TSSetFunction(), TSGetFunction()
2655 @*/
2656 PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
2657 {
2658   PetscErrorCode  ierr;
2659   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2660   int             nlhs = 2,nrhs = 9;
2661   mxArray         *plhs[2],*prhs[9];
2662   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
2663 
2664   PetscFunctionBegin;
2665   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2666   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2667 
2668   /* call Matlab function in ctx with arguments x and y */
2669 
2670   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
2671   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2672   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
2673   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
2674   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
2675   prhs[0] =  mxCreateDoubleScalar((double)ls);
2676   prhs[1] =  mxCreateDoubleScalar((double)time);
2677   prhs[2] =  mxCreateDoubleScalar((double)lx);
2678   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2679   prhs[4] =  mxCreateDoubleScalar((double)shift);
2680   prhs[5] =  mxCreateDoubleScalar((double)lA);
2681   prhs[6] =  mxCreateDoubleScalar((double)lB);
2682   prhs[7] =  mxCreateString(sctx->funcname);
2683   prhs[8] =  sctx->ctx;
2684   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
2685   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2686   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
2687   mxDestroyArray(prhs[0]);
2688   mxDestroyArray(prhs[1]);
2689   mxDestroyArray(prhs[2]);
2690   mxDestroyArray(prhs[3]);
2691   mxDestroyArray(prhs[4]);
2692   mxDestroyArray(prhs[5]);
2693   mxDestroyArray(prhs[6]);
2694   mxDestroyArray(prhs[7]);
2695   mxDestroyArray(plhs[0]);
2696   mxDestroyArray(plhs[1]);
2697   PetscFunctionReturn(0);
2698 }
2699 
2700 
2701 #undef __FUNCT__
2702 #define __FUNCT__ "TSSetJacobianMatlab"
2703 /*
2704    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
2705    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
2706 
2707    Logically Collective on TS
2708 
2709    Input Parameters:
2710 +  ts - the TS context
2711 .  A,B - Jacobian matrices
2712 .  func - function evaluation routine
2713 -  ctx - user context
2714 
2715    Calling sequence of func:
2716 $    flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
2717 
2718 
2719    Level: developer
2720 
2721 .keywords: TS, nonlinear, set, function
2722 
2723 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2724 */
2725 PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
2726 {
2727   PetscErrorCode    ierr;
2728   TSMatlabContext *sctx;
2729 
2730   PetscFunctionBegin;
2731   /* currently sctx is memory bleed */
2732   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2733   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2734   /*
2735      This should work, but it doesn't
2736   sctx->ctx = ctx;
2737   mexMakeArrayPersistent(sctx->ctx);
2738   */
2739   sctx->ctx = mxDuplicateArray(ctx);
2740   ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
2741   PetscFunctionReturn(0);
2742 }
2743 
2744 #undef __FUNCT__
2745 #define __FUNCT__ "TSMonitor_Matlab"
2746 /*
2747    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
2748 
2749    Collective on TS
2750 
2751 .seealso: TSSetFunction(), TSGetFunction()
2752 @*/
2753 PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
2754 {
2755   PetscErrorCode  ierr;
2756   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2757   int             nlhs = 1,nrhs = 6;
2758   mxArray         *plhs[1],*prhs[6];
2759   long long int   lx = 0,ls = 0;
2760 
2761   PetscFunctionBegin;
2762   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2763   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
2764 
2765   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
2766   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2767   prhs[0] =  mxCreateDoubleScalar((double)ls);
2768   prhs[1] =  mxCreateDoubleScalar((double)it);
2769   prhs[2] =  mxCreateDoubleScalar((double)time);
2770   prhs[3] =  mxCreateDoubleScalar((double)lx);
2771   prhs[4] =  mxCreateString(sctx->funcname);
2772   prhs[5] =  sctx->ctx;
2773   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
2774   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2775   mxDestroyArray(prhs[0]);
2776   mxDestroyArray(prhs[1]);
2777   mxDestroyArray(prhs[2]);
2778   mxDestroyArray(prhs[3]);
2779   mxDestroyArray(prhs[4]);
2780   mxDestroyArray(plhs[0]);
2781   PetscFunctionReturn(0);
2782 }
2783 
2784 
2785 #undef __FUNCT__
2786 #define __FUNCT__ "TSMonitorSetMatlab"
2787 /*
2788    TSMonitorSetMatlab - Sets the monitor function from Matlab
2789 
2790    Level: developer
2791 
2792 .keywords: TS, nonlinear, set, function
2793 
2794 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2795 */
2796 PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
2797 {
2798   PetscErrorCode    ierr;
2799   TSMatlabContext *sctx;
2800 
2801   PetscFunctionBegin;
2802   /* currently sctx is memory bleed */
2803   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2804   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2805   /*
2806      This should work, but it doesn't
2807   sctx->ctx = ctx;
2808   mexMakeArrayPersistent(sctx->ctx);
2809   */
2810   sctx->ctx = mxDuplicateArray(ctx);
2811   ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
2812   PetscFunctionReturn(0);
2813 }
2814 
2815 #endif
2816