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