xref: /petsc/src/ts/interface/ts.c (revision f3fe499b4cc4d64bf04aa4f5e4963dcc4eb56541)
1 #define PETSCTS_DLL
2 
3 #include "private/tsimpl.h"        /*I "petscts.h"  I*/
4 
5 /* Logging support */
6 PetscClassId PETSCTS_DLLEXPORT TS_CLASSID;
7 PetscLogEvent  TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "TSSetTypeFromOptions"
11 /*
12   TSSetTypeFromOptions - Sets the type of ts from user options.
13 
14   Collective on TS
15 
16   Input Parameter:
17 . ts - The ts
18 
19   Level: intermediate
20 
21 .keywords: TS, set, options, database, type
22 .seealso: TSSetFromOptions(), TSSetType()
23 */
24 static PetscErrorCode TSSetTypeFromOptions(TS ts)
25 {
26   PetscTruth     opt;
27   const char     *defaultType;
28   char           typeName[256];
29   PetscErrorCode ierr;
30 
31   PetscFunctionBegin;
32   if (((PetscObject)ts)->type_name) {
33     defaultType = ((PetscObject)ts)->type_name;
34   } else {
35     defaultType = TSEULER;
36   }
37 
38   if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
39   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
40   if (opt) {
41     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
42   } else {
43     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
44   }
45   PetscFunctionReturn(0);
46 }
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "TSSetFromOptions"
50 /*@
51    TSSetFromOptions - Sets various TS parameters from user options.
52 
53    Collective on TS
54 
55    Input Parameter:
56 .  ts - the TS context obtained from TSCreate()
57 
58    Options Database Keys:
59 +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP
60 .  -ts_max_steps maxsteps - maximum number of time-steps to take
61 .  -ts_max_time time - maximum time to compute to
62 .  -ts_dt dt - initial time step
63 .  -ts_monitor - print information at each timestep
64 -  -ts_monitor_draw - plot information at each timestep
65 
66    Level: beginner
67 
68 .keywords: TS, timestep, set, options, database
69 
70 .seealso: TSGetType()
71 @*/
72 PetscErrorCode PETSCTS_DLLEXPORT TSSetFromOptions(TS ts)
73 {
74   PetscReal               dt;
75   PetscTruth              opt,flg;
76   PetscErrorCode          ierr;
77   PetscViewerASCIIMonitor monviewer;
78   char                    monfilename[PETSC_MAX_PATH_LEN];
79 
80   PetscFunctionBegin;
81   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
82   ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr);
83 
84     /* Handle generic TS options */
85     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
86     ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
87     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
88     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
89     if (opt) {
90       ts->initial_time_step = ts->time_step = dt;
91     }
92 
93     /* Monitor options */
94     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
95     if (flg) {
96       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,monfilename,((PetscObject)ts)->tablevel,&monviewer);CHKERRQ(ierr);
97       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
98     }
99     opt  = PETSC_FALSE;
100     ierr = PetscOptionsTruth("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
101     if (opt) {
102       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
103     }
104     opt  = PETSC_FALSE;
105     ierr = PetscOptionsTruth("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
106     if (opt) {
107       ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
108     }
109 
110     /* Handle TS type options */
111     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
112 
113     /* Handle specific TS options */
114     if (ts->ops->setfromoptions) {
115       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
116     }
117 
118     /* process any options handlers added with PetscObjectAddOptionsHandler() */
119     ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
120   ierr = PetscOptionsEnd();CHKERRQ(ierr);
121 
122   /* Handle subobject options */
123   switch(ts->problem_type) {
124     /* Should check for implicit/explicit */
125   case TS_LINEAR:
126     if (ts->ksp) {
127       ierr = KSPSetOperators(ts->ksp,ts->Arhs,ts->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
128       ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr);
129     }
130     break;
131   case TS_NONLINEAR:
132     if (ts->snes) {
133       /* this is a bit of a hack, but it gets the matrix information into SNES earlier
134          so that SNES and KSP have more information to pick reasonable defaults
135          before they allow users to set options
136        * If ts->A has been set at this point, we are probably using the implicit form
137          and Arhs will never be used. */
138       ierr = SNESSetJacobian(ts->snes,ts->A?ts->A:ts->Arhs,ts->B,0,ts);CHKERRQ(ierr);
139       ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
140     }
141     break;
142   default:
143     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type);
144   }
145 
146   PetscFunctionReturn(0);
147 }
148 
149 #undef  __FUNCT__
150 #define __FUNCT__ "TSViewFromOptions"
151 /*@
152   TSViewFromOptions - This function visualizes the ts based upon user options.
153 
154   Collective on TS
155 
156   Input Parameter:
157 . ts - The ts
158 
159   Level: intermediate
160 
161 .keywords: TS, view, options, database
162 .seealso: TSSetFromOptions(), TSView()
163 @*/
164 PetscErrorCode PETSCTS_DLLEXPORT TSViewFromOptions(TS ts,const char title[])
165 {
166   PetscViewer    viewer;
167   PetscDraw      draw;
168   PetscTruth     opt = PETSC_FALSE;
169   char           fileName[PETSC_MAX_PATH_LEN];
170   PetscErrorCode ierr;
171 
172   PetscFunctionBegin;
173   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr);
174   if (opt && !PetscPreLoadingOn) {
175     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr);
176     ierr = TSView(ts, viewer);CHKERRQ(ierr);
177     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
178   }
179   opt = PETSC_FALSE;
180   ierr = PetscOptionsGetTruth(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr);
181   if (opt) {
182     ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
183     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
184     if (title) {
185       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
186     } else {
187       ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr);
188       ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr);
189     }
190     ierr = TSView(ts, viewer);CHKERRQ(ierr);
191     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
192     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
193     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
194   }
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "TSComputeRHSJacobian"
200 /*@
201    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
202       set with TSSetRHSJacobian().
203 
204    Collective on TS and Vec
205 
206    Input Parameters:
207 +  ts - the TS context
208 .  t - current timestep
209 -  x - input vector
210 
211    Output Parameters:
212 +  A - Jacobian matrix
213 .  B - optional preconditioning matrix
214 -  flag - flag indicating matrix structure
215 
216    Notes:
217    Most users should not need to explicitly call this routine, as it
218    is used internally within the nonlinear solvers.
219 
220    See KSPSetOperators() for important information about setting the
221    flag parameter.
222 
223    TSComputeJacobian() is valid only for TS_NONLINEAR
224 
225    Level: developer
226 
227 .keywords: SNES, compute, Jacobian, matrix
228 
229 .seealso:  TSSetRHSJacobian(), KSPSetOperators()
230 @*/
231 PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
232 {
233   PetscErrorCode ierr;
234 
235   PetscFunctionBegin;
236   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
237   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
238   PetscCheckSameComm(ts,1,X,3);
239   if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
240   if (ts->ops->rhsjacobian) {
241     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
242     *flg = DIFFERENT_NONZERO_PATTERN;
243     PetscStackPush("TS user Jacobian function");
244     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
245     PetscStackPop;
246     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
247     /* make sure user returned a correct Jacobian and preconditioner */
248     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
249     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
250   } else {
251     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
253     if (*A != *B) {
254       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256     }
257   }
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "TSComputeRHSFunction"
263 /*@
264    TSComputeRHSFunction - Evaluates the right-hand-side function.
265 
266    Collective on TS and Vec
267 
268    Input Parameters:
269 +  ts - the TS context
270 .  t - current time
271 -  x - state vector
272 
273    Output Parameter:
274 .  y - right hand side
275 
276    Note:
277    Most users should not need to explicitly call this routine, as it
278    is used internally within the nonlinear solvers.
279 
280    If the user did not provide a function but merely a matrix,
281    this routine applies the matrix.
282 
283    Level: developer
284 
285 .keywords: TS, compute
286 
287 .seealso: TSSetRHSFunction(), TSComputeIFunction()
288 @*/
289 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
290 {
291   PetscErrorCode ierr;
292 
293   PetscFunctionBegin;
294   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
295   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
296   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
297 
298   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
299   if (ts->ops->rhsfunction) {
300     PetscStackPush("TS user right-hand-side function");
301     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
302     PetscStackPop;
303   } else {
304     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
305       MatStructure flg;
306       PetscStackPush("TS user right-hand-side matrix function");
307       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
308       PetscStackPop;
309     }
310     ierr = MatMult(ts->Arhs,x,y);CHKERRQ(ierr);
311   }
312 
313   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
314 
315   PetscFunctionReturn(0);
316 }
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "TSComputeIFunction"
320 /*@
321    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
322 
323    Collective on TS and Vec
324 
325    Input Parameters:
326 +  ts - the TS context
327 .  t - current time
328 .  X - state vector
329 -  Xdot - time derivative of state vector
330 
331    Output Parameter:
332 .  Y - right hand side
333 
334    Note:
335    Most users should not need to explicitly call this routine, as it
336    is used internally within the nonlinear solvers.
337 
338    If the user did did not write their equations in implicit form, this
339    function recasts them in implicit form.
340 
341    Level: developer
342 
343 .keywords: TS, compute
344 
345 .seealso: TSSetIFunction(), TSComputeRHSFunction()
346 @*/
347 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y)
348 {
349   PetscErrorCode ierr;
350 
351   PetscFunctionBegin;
352   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
353   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
354   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
355   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
356 
357   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
358   if (ts->ops->ifunction) {
359     PetscStackPush("TS user implicit function");
360     ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr);
361     PetscStackPop;
362   } else {
363     if (ts->ops->rhsfunction) {
364       PetscStackPush("TS user right-hand-side function");
365       ierr = (*ts->ops->rhsfunction)(ts,t,X,Y,ts->funP);CHKERRQ(ierr);
366       PetscStackPop;
367     } else {
368       if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
369         MatStructure flg;
370         /* Note: flg is not being used.
371            For it to be useful, we'd have to cache it and then apply it in TSComputeIJacobian.
372         */
373         PetscStackPush("TS user right-hand-side matrix function");
374         ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
375         PetscStackPop;
376       }
377       ierr = MatMult(ts->Arhs,X,Y);CHKERRQ(ierr);
378     }
379 
380     /* Convert to implicit form: F(X,Xdot) = Alhs * Xdot - Frhs(X) */
381     if (ts->Alhs) {
382       if (ts->ops->lhsmatrix) {
383         MatStructure flg;
384         PetscStackPush("TS user left-hand-side matrix function");
385         ierr = (*ts->ops->lhsmatrix)(ts,t,&ts->Alhs,PETSC_NULL,&flg,ts->jacP);CHKERRQ(ierr);
386         PetscStackPop;
387       }
388       ierr = VecScale(Y,-1.);CHKERRQ(ierr);
389       ierr = MatMultAdd(ts->Alhs,Xdot,Y,Y);CHKERRQ(ierr);
390     } else {
391       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
392     }
393   }
394   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNCT__
399 #define __FUNCT__ "TSComputeIJacobian"
400 /*@
401    TSComputeIJacobian - Evaluates the Jacobian of the DAE
402 
403    Collective on TS and Vec
404 
405    Input
406       Input Parameters:
407 +  ts - the TS context
408 .  t - current timestep
409 .  X - state vector
410 .  Xdot - time derivative of state vector
411 -  shift - shift to apply, see note below
412 
413    Output Parameters:
414 +  A - Jacobian matrix
415 .  B - optional preconditioning matrix
416 -  flag - flag indicating matrix structure
417 
418    Notes:
419    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
420 
421    dF/dX + shift*dF/dXdot
422 
423    Most users should not need to explicitly call this routine, as it
424    is used internally within the nonlinear solvers.
425 
426    TSComputeIJacobian() is valid only for TS_NONLINEAR
427 
428    Level: developer
429 
430 .keywords: TS, compute, Jacobian, matrix
431 
432 .seealso:  TSSetIJacobian()
433 @*/
434 PetscErrorCode PETSCTS_DLLEXPORT TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg)
435 {
436   PetscErrorCode ierr;
437 
438   PetscFunctionBegin;
439   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
440   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
441   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
442   PetscValidPointer(A,6);
443   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
444   PetscValidPointer(B,7);
445   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
446   PetscValidPointer(flg,8);
447 
448   *flg = SAME_NONZERO_PATTERN;  /* In case it we're solving a linear problem in which case it wouldn't get initialized below. */
449   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
450   if (ts->ops->ijacobian) {
451     PetscStackPush("TS user implicit Jacobian");
452     ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr);
453     PetscStackPop;
454   } else {
455     if (ts->ops->rhsjacobian) {
456       PetscStackPush("TS user right-hand-side Jacobian");
457       ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
458       PetscStackPop;
459     } else {
460       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
461       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
462       if (*A != *B) {
463         ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
464         ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
465       }
466     }
467 
468     /* Convert to implicit form */
469     /* inefficient because these operations will normally traverse all matrix elements twice */
470     ierr = MatScale(*A,-1);CHKERRQ(ierr);
471     if (ts->Alhs) {
472       ierr = MatAXPY(*A,shift,ts->Alhs,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
473     } else {
474       ierr = MatShift(*A,shift);CHKERRQ(ierr);
475     }
476     if (*A != *B) {
477       ierr = MatScale(*B,-1);CHKERRQ(ierr);
478       ierr = MatShift(*B,shift);CHKERRQ(ierr);
479     }
480   }
481   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
482   PetscFunctionReturn(0);
483 }
484 
485 #undef __FUNCT__
486 #define __FUNCT__ "TSSetRHSFunction"
487 /*@C
488     TSSetRHSFunction - Sets the routine for evaluating the function,
489     F(t,u), where U_t = F(t,u).
490 
491     Logically Collective on TS
492 
493     Input Parameters:
494 +   ts - the TS context obtained from TSCreate()
495 .   f - routine for evaluating the right-hand-side function
496 -   ctx - [optional] user-defined context for private data for the
497           function evaluation routine (may be PETSC_NULL)
498 
499     Calling sequence of func:
500 $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
501 
502 +   t - current timestep
503 .   u - input vector
504 .   F - function vector
505 -   ctx - [optional] user-defined function context
506 
507     Important:
508     The user MUST call either this routine or TSSetMatrices().
509 
510     Level: beginner
511 
512 .keywords: TS, timestep, set, right-hand-side, function
513 
514 .seealso: TSSetMatrices()
515 @*/
516 PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
517 {
518   PetscFunctionBegin;
519 
520   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
521   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
522   ts->ops->rhsfunction = f;
523   ts->funP             = ctx;
524   PetscFunctionReturn(0);
525 }
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "TSSetMatrices"
529 /*@C
530    TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs,
531    where Alhs(t) U_t = Arhs(t) U.
532 
533    Logically Collective on TS
534 
535    Input Parameters:
536 +  ts   - the TS context obtained from TSCreate()
537 .  Arhs - matrix
538 .  frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
539           if Arhs is not a function of t.
540 .  Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix.
541 .  flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
542           if Alhs is not a function of t.
543 .  flag - flag indicating information about the matrix structure of Arhs and Alhs.
544           The available options are
545             SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs
546             DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs
547 -  ctx  - [optional] user-defined context for private data for the
548           matrix evaluation routine (may be PETSC_NULL)
549 
550    Calling sequence of func:
551 $     func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
552 
553 +  t - current timestep
554 .  A - matrix A, where U_t = A(t) U
555 .  B - preconditioner matrix, usually the same as A
556 .  flag - flag indicating information about the preconditioner matrix
557           structure (same as flag in KSPSetOperators())
558 -  ctx - [optional] user-defined context for matrix evaluation routine
559 
560    Notes:
561    The routine func() takes Mat* as the matrix arguments rather than Mat.
562    This allows the matrix evaluation routine to replace Arhs or Alhs with a
563    completely new new matrix structure (not just different matrix elements)
564    when appropriate, for instance, if the nonzero structure is changing
565    throughout the global iterations.
566 
567    Important:
568    The user MUST call either this routine or TSSetRHSFunction().
569 
570    Level: beginner
571 
572 .keywords: TS, timestep, set, matrix
573 
574 .seealso: TSSetRHSFunction()
575 @*/
576 PetscErrorCode PETSCTS_DLLEXPORT TSSetMatrices(TS ts,Mat Arhs,PetscErrorCode (*frhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat Alhs,PetscErrorCode (*flhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure flag,void *ctx)
577 {
578   PetscFunctionBegin;
579   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
580   if (Arhs){
581     PetscValidHeaderSpecific(Arhs,MAT_CLASSID,2);
582     PetscCheckSameComm(ts,1,Arhs,2);
583     ts->Arhs           = Arhs;
584     ts->ops->rhsmatrix = frhs;
585   }
586   if (Alhs){
587     PetscValidHeaderSpecific(Alhs,MAT_CLASSID,4);
588     PetscCheckSameComm(ts,1,Alhs,4);
589     ts->Alhs           = Alhs;
590     ts->ops->lhsmatrix = flhs;
591   }
592 
593   ts->jacP           = ctx;
594   ts->matflg         = flag;
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "TSGetMatrices"
600 /*@C
601    TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep,
602    where Alhs(t) U_t = Arhs(t) U.
603 
604    Not Collective, but parallel objects are returned if TS is parallel
605 
606    Input Parameter:
607 .  ts  - The TS context obtained from TSCreate()
608 
609    Output Parameters:
610 +  Arhs - The right-hand side matrix
611 .  Alhs - The left-hand side matrix
612 -  ctx - User-defined context for matrix evaluation routine
613 
614    Notes: You can pass in PETSC_NULL for any return argument you do not need.
615 
616    Level: intermediate
617 
618 .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
619 
620 .keywords: TS, timestep, get, matrix
621 
622 @*/
623 PetscErrorCode PETSCTS_DLLEXPORT TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx)
624 {
625   PetscFunctionBegin;
626   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
627   if (Arhs) *Arhs = ts->Arhs;
628   if (Alhs) *Alhs = ts->Alhs;
629   if (ctx)  *ctx = ts->jacP;
630   PetscFunctionReturn(0);
631 }
632 
633 #undef __FUNCT__
634 #define __FUNCT__ "TSSetRHSJacobian"
635 /*@C
636    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
637    where U_t = F(U,t), as well as the location to store the matrix.
638    Use TSSetMatrices() for linear problems.
639 
640    Logically Collective on TS
641 
642    Input Parameters:
643 +  ts  - the TS context obtained from TSCreate()
644 .  A   - Jacobian matrix
645 .  B   - preconditioner matrix (usually same as A)
646 .  f   - the Jacobian evaluation routine
647 -  ctx - [optional] user-defined context for private data for the
648          Jacobian evaluation routine (may be PETSC_NULL)
649 
650    Calling sequence of func:
651 $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
652 
653 +  t - current timestep
654 .  u - input vector
655 .  A - matrix A, where U_t = A(t)u
656 .  B - preconditioner matrix, usually the same as A
657 .  flag - flag indicating information about the preconditioner matrix
658           structure (same as flag in KSPSetOperators())
659 -  ctx - [optional] user-defined context for matrix evaluation routine
660 
661    Notes:
662    See KSPSetOperators() for important information about setting the flag
663    output parameter in the routine func().  Be sure to read this information!
664 
665    The routine func() takes Mat * as the matrix arguments rather than Mat.
666    This allows the matrix evaluation routine to replace A and/or B with a
667    completely new matrix structure (not just different matrix elements)
668    when appropriate, for instance, if the nonzero structure is changing
669    throughout the global iterations.
670 
671    Level: beginner
672 
673 .keywords: TS, timestep, set, right-hand-side, Jacobian
674 
675 .seealso: TSDefaultComputeJacobianColor(),
676           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices()
677 
678 @*/
679 PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
680 {
681   PetscFunctionBegin;
682   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
683   PetscValidHeaderSpecific(A,MAT_CLASSID,2);
684   PetscValidHeaderSpecific(B,MAT_CLASSID,3);
685   PetscCheckSameComm(ts,1,A,2);
686   PetscCheckSameComm(ts,1,B,3);
687   if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()");
688 
689   ts->ops->rhsjacobian = f;
690   ts->jacP             = ctx;
691   ts->Arhs             = A;
692   ts->B                = B;
693   PetscFunctionReturn(0);
694 }
695 
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "TSSetIFunction"
699 /*@C
700    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
701 
702    Logically Collective on TS
703 
704    Input Parameters:
705 +  ts  - the TS context obtained from TSCreate()
706 .  f   - the function evaluation routine
707 -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
708 
709    Calling sequence of f:
710 $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
711 
712 +  t   - time at step/stage being solved
713 .  u   - state vector
714 .  u_t - time derivative of state vector
715 .  F   - function vector
716 -  ctx - [optional] user-defined context for matrix evaluation routine
717 
718    Important:
719    The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices().  This routine must be used when not solving an ODE.
720 
721    Level: beginner
722 
723 .keywords: TS, timestep, set, DAE, Jacobian
724 
725 .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian()
726 @*/
727 PetscErrorCode PETSCTS_DLLEXPORT TSSetIFunction(TS ts,TSIFunction f,void *ctx)
728 {
729 
730   PetscFunctionBegin;
731   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
732   ts->ops->ifunction = f;
733   ts->funP           = ctx;
734   PetscFunctionReturn(0);
735 }
736 
737 #undef __FUNCT__
738 #define __FUNCT__ "TSSetIJacobian"
739 /*@C
740    TSSetIJacobian - Set the function to compute the Jacobian of
741    G(U) = F(t,U,U0+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
742 
743    Logically Collective on TS
744 
745    Input Parameters:
746 +  ts  - the TS context obtained from TSCreate()
747 .  A   - Jacobian matrix
748 .  B   - preconditioning matrix for A (may be same as A)
749 .  f   - the Jacobian evaluation routine
750 -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
751 
752    Calling sequence of f:
753 $  f(TS ts,PetscReal t,Vec u,Vec u_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
754 
755 +  t    - time at step/stage being solved
756 .  u    - state vector
757 .  u_t  - time derivative of state vector
758 .  a    - shift
759 .  A    - Jacobian of G(U) = F(t,U,U0+a*U), equivalent to dF/dU + a*dF/dU_t
760 .  B    - preconditioning matrix for A, may be same as A
761 .  flag - flag indicating information about the preconditioner matrix
762           structure (same as flag in KSPSetOperators())
763 -  ctx  - [optional] user-defined context for matrix evaluation routine
764 
765    Notes:
766    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
767 
768    Level: beginner
769 
770 .keywords: TS, timestep, DAE, Jacobian
771 
772 .seealso: TSSetIFunction(), TSSetRHSJacobian()
773 
774 @*/
775 PetscErrorCode PETSCTS_DLLEXPORT TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
776 {
777   PetscErrorCode ierr;
778 
779   PetscFunctionBegin;
780   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
781   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
782   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
783   if (A) PetscCheckSameComm(ts,1,A,2);
784   if (B) PetscCheckSameComm(ts,1,B,3);
785   if (f)   ts->ops->ijacobian = f;
786   if (ctx) ts->jacP             = ctx;
787   if (A) {
788     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
789     if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);}
790     ts->A = A;
791   }
792 #if 0
793   /* The sane and consistent alternative */
794   if (B) {
795     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
796     if (ts->B) {ierr = MatDestroy(ts->B);CHKERRQ(ierr);}
797     ts->B = B;
798   }
799 #else
800   /* Don't reference B because TSDestroy() doesn't destroy it.  These ownership semantics are awkward and inconsistent. */
801   if (B) ts->B = B;
802 #endif
803   PetscFunctionReturn(0);
804 }
805 
806 #undef __FUNCT__
807 #define __FUNCT__ "TSView"
808 /*@C
809     TSView - Prints the TS data structure.
810 
811     Collective on TS
812 
813     Input Parameters:
814 +   ts - the TS context obtained from TSCreate()
815 -   viewer - visualization context
816 
817     Options Database Key:
818 .   -ts_view - calls TSView() at end of TSStep()
819 
820     Notes:
821     The available visualization contexts include
822 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
823 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
824          output where only the first processor opens
825          the file.  All other processors send their
826          data to the first processor to print.
827 
828     The user can open an alternative visualization context with
829     PetscViewerASCIIOpen() - output to a specified file.
830 
831     Level: beginner
832 
833 .keywords: TS, timestep, view
834 
835 .seealso: PetscViewerASCIIOpen()
836 @*/
837 PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer)
838 {
839   PetscErrorCode ierr;
840   const TSType   type;
841   PetscTruth     iascii,isstring;
842 
843   PetscFunctionBegin;
844   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
845   if (!viewer) {
846     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
847   }
848   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
849   PetscCheckSameComm(ts,1,viewer,2);
850 
851   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
852   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
853   if (iascii) {
854     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
855     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
856     if (type) {
857       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
858     } else {
859       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
860     }
861     if (ts->ops->view) {
862       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
863       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
864       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
865     }
866     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
867     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
868     if (ts->problem_type == TS_NONLINEAR) {
869       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
870     }
871     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
872   } else if (isstring) {
873     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
874     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
875   }
876   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
877   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
878   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
879   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
880   PetscFunctionReturn(0);
881 }
882 
883 
884 #undef __FUNCT__
885 #define __FUNCT__ "TSSetApplicationContext"
886 /*@C
887    TSSetApplicationContext - Sets an optional user-defined context for
888    the timesteppers.
889 
890    Logically Collective on TS
891 
892    Input Parameters:
893 +  ts - the TS context obtained from TSCreate()
894 -  usrP - optional user context
895 
896    Level: intermediate
897 
898 .keywords: TS, timestep, set, application, context
899 
900 .seealso: TSGetApplicationContext()
901 @*/
902 PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP)
903 {
904   PetscFunctionBegin;
905   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
906   ts->user = usrP;
907   PetscFunctionReturn(0);
908 }
909 
910 #undef __FUNCT__
911 #define __FUNCT__ "TSGetApplicationContext"
912 /*@C
913     TSGetApplicationContext - Gets the user-defined context for the
914     timestepper.
915 
916     Not Collective
917 
918     Input Parameter:
919 .   ts - the TS context obtained from TSCreate()
920 
921     Output Parameter:
922 .   usrP - user context
923 
924     Level: intermediate
925 
926 .keywords: TS, timestep, get, application, context
927 
928 .seealso: TSSetApplicationContext()
929 @*/
930 PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP)
931 {
932   PetscFunctionBegin;
933   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
934   *usrP = ts->user;
935   PetscFunctionReturn(0);
936 }
937 
938 #undef __FUNCT__
939 #define __FUNCT__ "TSGetTimeStepNumber"
940 /*@
941    TSGetTimeStepNumber - Gets the current number of timesteps.
942 
943    Not Collective
944 
945    Input Parameter:
946 .  ts - the TS context obtained from TSCreate()
947 
948    Output Parameter:
949 .  iter - number steps so far
950 
951    Level: intermediate
952 
953 .keywords: TS, timestep, get, iteration, number
954 @*/
955 PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter)
956 {
957   PetscFunctionBegin;
958   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
959   PetscValidIntPointer(iter,2);
960   *iter = ts->steps;
961   PetscFunctionReturn(0);
962 }
963 
964 #undef __FUNCT__
965 #define __FUNCT__ "TSSetInitialTimeStep"
966 /*@
967    TSSetInitialTimeStep - Sets the initial timestep to be used,
968    as well as the initial time.
969 
970    Logically Collective on TS
971 
972    Input Parameters:
973 +  ts - the TS context obtained from TSCreate()
974 .  initial_time - the initial time
975 -  time_step - the size of the timestep
976 
977    Level: intermediate
978 
979 .seealso: TSSetTimeStep(), TSGetTimeStep()
980 
981 .keywords: TS, set, initial, timestep
982 @*/
983 PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
984 {
985   PetscFunctionBegin;
986   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
987   ts->time_step         = time_step;
988   ts->initial_time_step = time_step;
989   ts->ptime             = initial_time;
990   PetscFunctionReturn(0);
991 }
992 
993 #undef __FUNCT__
994 #define __FUNCT__ "TSSetTimeStep"
995 /*@
996    TSSetTimeStep - Allows one to reset the timestep at any time,
997    useful for simple pseudo-timestepping codes.
998 
999    Logically Collective on TS
1000 
1001    Input Parameters:
1002 +  ts - the TS context obtained from TSCreate()
1003 -  time_step - the size of the timestep
1004 
1005    Level: intermediate
1006 
1007 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1008 
1009 .keywords: TS, set, timestep
1010 @*/
1011 PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step)
1012 {
1013   PetscFunctionBegin;
1014   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1015   PetscValidLogicalCollectiveReal(ts,time_step,2);
1016   ts->time_step = time_step;
1017   PetscFunctionReturn(0);
1018 }
1019 
1020 #undef __FUNCT__
1021 #define __FUNCT__ "TSGetTimeStep"
1022 /*@
1023    TSGetTimeStep - Gets the current timestep size.
1024 
1025    Not Collective
1026 
1027    Input Parameter:
1028 .  ts - the TS context obtained from TSCreate()
1029 
1030    Output Parameter:
1031 .  dt - the current timestep size
1032 
1033    Level: intermediate
1034 
1035 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1036 
1037 .keywords: TS, get, timestep
1038 @*/
1039 PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt)
1040 {
1041   PetscFunctionBegin;
1042   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1043   PetscValidDoublePointer(dt,2);
1044   *dt = ts->time_step;
1045   PetscFunctionReturn(0);
1046 }
1047 
1048 #undef __FUNCT__
1049 #define __FUNCT__ "TSGetSolution"
1050 /*@
1051    TSGetSolution - Returns the solution at the present timestep. It
1052    is valid to call this routine inside the function that you are evaluating
1053    in order to move to the new timestep. This vector not changed until
1054    the solution at the next timestep has been calculated.
1055 
1056    Not Collective, but Vec returned is parallel if TS is parallel
1057 
1058    Input Parameter:
1059 .  ts - the TS context obtained from TSCreate()
1060 
1061    Output Parameter:
1062 .  v - the vector containing the solution
1063 
1064    Level: intermediate
1065 
1066 .seealso: TSGetTimeStep()
1067 
1068 .keywords: TS, timestep, get, solution
1069 @*/
1070 PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v)
1071 {
1072   PetscFunctionBegin;
1073   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1074   PetscValidPointer(v,2);
1075   *v = ts->vec_sol_always;
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 /* ----- Routines to initialize and destroy a timestepper ---- */
1080 #undef __FUNCT__
1081 #define __FUNCT__ "TSSetProblemType"
1082 /*@
1083   TSSetProblemType - Sets the type of problem to be solved.
1084 
1085   Not collective
1086 
1087   Input Parameters:
1088 + ts   - The TS
1089 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1090 .vb
1091          U_t = A U
1092          U_t = A(t) U
1093          U_t = F(t,U)
1094 .ve
1095 
1096    Level: beginner
1097 
1098 .keywords: TS, problem type
1099 .seealso: TSSetUp(), TSProblemType, TS
1100 @*/
1101 PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type)
1102 {
1103   PetscFunctionBegin;
1104   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1105   ts->problem_type = type;
1106   PetscFunctionReturn(0);
1107 }
1108 
1109 #undef __FUNCT__
1110 #define __FUNCT__ "TSGetProblemType"
1111 /*@C
1112   TSGetProblemType - Gets the type of problem to be solved.
1113 
1114   Not collective
1115 
1116   Input Parameter:
1117 . ts   - The TS
1118 
1119   Output Parameter:
1120 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1121 .vb
1122          U_t = A U
1123          U_t = A(t) U
1124          U_t = F(t,U)
1125 .ve
1126 
1127    Level: beginner
1128 
1129 .keywords: TS, problem type
1130 .seealso: TSSetUp(), TSProblemType, TS
1131 @*/
1132 PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type)
1133 {
1134   PetscFunctionBegin;
1135   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1136   PetscValidIntPointer(type,2);
1137   *type = ts->problem_type;
1138   PetscFunctionReturn(0);
1139 }
1140 
1141 #undef __FUNCT__
1142 #define __FUNCT__ "TSSetUp"
1143 /*@
1144    TSSetUp - Sets up the internal data structures for the later use
1145    of a timestepper.
1146 
1147    Collective on TS
1148 
1149    Input Parameter:
1150 .  ts - the TS context obtained from TSCreate()
1151 
1152    Notes:
1153    For basic use of the TS solvers the user need not explicitly call
1154    TSSetUp(), since these actions will automatically occur during
1155    the call to TSStep().  However, if one wishes to control this
1156    phase separately, TSSetUp() should be called after TSCreate()
1157    and optional routines of the form TSSetXXX(), but before TSStep().
1158 
1159    Level: advanced
1160 
1161 .keywords: TS, timestep, setup
1162 
1163 .seealso: TSCreate(), TSStep(), TSDestroy()
1164 @*/
1165 PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts)
1166 {
1167   PetscErrorCode ierr;
1168 
1169   PetscFunctionBegin;
1170   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1171   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1172   if (!((PetscObject)ts)->type_name) {
1173     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1174   }
1175   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1176   ts->setupcalled = 1;
1177   PetscFunctionReturn(0);
1178 }
1179 
1180 #undef __FUNCT__
1181 #define __FUNCT__ "TSDestroy"
1182 /*@
1183    TSDestroy - Destroys the timestepper context that was created
1184    with TSCreate().
1185 
1186    Collective on TS
1187 
1188    Input Parameter:
1189 .  ts - the TS context obtained from TSCreate()
1190 
1191    Level: beginner
1192 
1193 .keywords: TS, timestepper, destroy
1194 
1195 .seealso: TSCreate(), TSSetUp(), TSSolve()
1196 @*/
1197 PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts)
1198 {
1199   PetscErrorCode ierr;
1200 
1201   PetscFunctionBegin;
1202   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1203   if (--((PetscObject)ts)->refct > 0) PetscFunctionReturn(0);
1204 
1205   /* if memory was published with AMS then destroy it */
1206   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
1207 
1208   if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);}
1209   if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);}
1210   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
1211   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
1212   if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);}
1213   ierr = TSMonitorCancel(ts);CHKERRQ(ierr);
1214   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1215   PetscFunctionReturn(0);
1216 }
1217 
1218 #undef __FUNCT__
1219 #define __FUNCT__ "TSGetSNES"
1220 /*@
1221    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1222    a TS (timestepper) context. Valid only for nonlinear problems.
1223 
1224    Not Collective, but SNES is parallel if TS is parallel
1225 
1226    Input Parameter:
1227 .  ts - the TS context obtained from TSCreate()
1228 
1229    Output Parameter:
1230 .  snes - the nonlinear solver context
1231 
1232    Notes:
1233    The user can then directly manipulate the SNES context to set various
1234    options, etc.  Likewise, the user can then extract and manipulate the
1235    KSP, KSP, and PC contexts as well.
1236 
1237    TSGetSNES() does not work for integrators that do not use SNES; in
1238    this case TSGetSNES() returns PETSC_NULL in snes.
1239 
1240    Level: beginner
1241 
1242 .keywords: timestep, get, SNES
1243 @*/
1244 PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes)
1245 {
1246   PetscFunctionBegin;
1247   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1248   PetscValidPointer(snes,2);
1249   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first");
1250   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
1251   *snes = ts->snes;
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 #undef __FUNCT__
1256 #define __FUNCT__ "TSGetKSP"
1257 /*@
1258    TSGetKSP - Returns the KSP (linear solver) associated with
1259    a TS (timestepper) context.
1260 
1261    Not Collective, but KSP is parallel if TS is parallel
1262 
1263    Input Parameter:
1264 .  ts - the TS context obtained from TSCreate()
1265 
1266    Output Parameter:
1267 .  ksp - the nonlinear solver context
1268 
1269    Notes:
1270    The user can then directly manipulate the KSP context to set various
1271    options, etc.  Likewise, the user can then extract and manipulate the
1272    KSP and PC contexts as well.
1273 
1274    TSGetKSP() does not work for integrators that do not use KSP;
1275    in this case TSGetKSP() returns PETSC_NULL in ksp.
1276 
1277    Level: beginner
1278 
1279 .keywords: timestep, get, KSP
1280 @*/
1281 PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp)
1282 {
1283   PetscFunctionBegin;
1284   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1285   PetscValidPointer(ksp,2);
1286   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1287   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1288   *ksp = ts->ksp;
1289   PetscFunctionReturn(0);
1290 }
1291 
1292 /* ----------- Routines to set solver parameters ---------- */
1293 
1294 #undef __FUNCT__
1295 #define __FUNCT__ "TSGetDuration"
1296 /*@
1297    TSGetDuration - Gets the maximum number of timesteps to use and
1298    maximum time for iteration.
1299 
1300    Not Collective
1301 
1302    Input Parameters:
1303 +  ts       - the TS context obtained from TSCreate()
1304 .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1305 -  maxtime  - final time to iterate to, or PETSC_NULL
1306 
1307    Level: intermediate
1308 
1309 .keywords: TS, timestep, get, maximum, iterations, time
1310 @*/
1311 PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1312 {
1313   PetscFunctionBegin;
1314   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1315   if (maxsteps) {
1316     PetscValidIntPointer(maxsteps,2);
1317     *maxsteps = ts->max_steps;
1318   }
1319   if (maxtime ) {
1320     PetscValidScalarPointer(maxtime,3);
1321     *maxtime  = ts->max_time;
1322   }
1323   PetscFunctionReturn(0);
1324 }
1325 
1326 #undef __FUNCT__
1327 #define __FUNCT__ "TSSetDuration"
1328 /*@
1329    TSSetDuration - Sets the maximum number of timesteps to use and
1330    maximum time for iteration.
1331 
1332    Logically Collective on TS
1333 
1334    Input Parameters:
1335 +  ts - the TS context obtained from TSCreate()
1336 .  maxsteps - maximum number of iterations to use
1337 -  maxtime - final time to iterate to
1338 
1339    Options Database Keys:
1340 .  -ts_max_steps <maxsteps> - Sets maxsteps
1341 .  -ts_max_time <maxtime> - Sets maxtime
1342 
1343    Notes:
1344    The default maximum number of iterations is 5000. Default time is 5.0
1345 
1346    Level: intermediate
1347 
1348 .keywords: TS, timestep, set, maximum, iterations
1349 @*/
1350 PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1351 {
1352   PetscFunctionBegin;
1353   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1354   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1355   PetscValidLogicalCollectiveReal(ts,maxtime,2);
1356   ts->max_steps = maxsteps;
1357   ts->max_time  = maxtime;
1358   PetscFunctionReturn(0);
1359 }
1360 
1361 #undef __FUNCT__
1362 #define __FUNCT__ "TSSetSolution"
1363 /*@
1364    TSSetSolution - Sets the initial solution vector
1365    for use by the TS routines.
1366 
1367    Logically Collective on TS and Vec
1368 
1369    Input Parameters:
1370 +  ts - the TS context obtained from TSCreate()
1371 -  x - the solution vector
1372 
1373    Level: beginner
1374 
1375 .keywords: TS, timestep, set, solution, initial conditions
1376 @*/
1377 PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x)
1378 {
1379   PetscFunctionBegin;
1380   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1381   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1382   ts->vec_sol        = ts->vec_sol_always = x;
1383   PetscFunctionReturn(0);
1384 }
1385 
1386 #undef __FUNCT__
1387 #define __FUNCT__ "TSSetPreStep"
1388 /*@C
1389   TSSetPreStep - Sets the general-purpose function
1390   called once at the beginning of each time step.
1391 
1392   Logically Collective on TS
1393 
1394   Input Parameters:
1395 + ts   - The TS context obtained from TSCreate()
1396 - func - The function
1397 
1398   Calling sequence of func:
1399 . func (TS ts);
1400 
1401   Level: intermediate
1402 
1403 .keywords: TS, timestep
1404 @*/
1405 PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1406 {
1407   PetscFunctionBegin;
1408   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1409   ts->ops->prestep = func;
1410   PetscFunctionReturn(0);
1411 }
1412 
1413 #undef __FUNCT__
1414 #define __FUNCT__ "TSPreStep"
1415 /*@C
1416   TSPreStep - Runs the user-defined pre-step function.
1417 
1418   Collective on TS
1419 
1420   Input Parameters:
1421 . ts   - The TS context obtained from TSCreate()
1422 
1423   Notes:
1424   TSPreStep() is typically used within time stepping implementations,
1425   so most users would not generally call this routine themselves.
1426 
1427   Level: developer
1428 
1429 .keywords: TS, timestep
1430 @*/
1431 PetscErrorCode PETSCTS_DLLEXPORT TSPreStep(TS ts)
1432 {
1433   PetscErrorCode ierr;
1434 
1435   PetscFunctionBegin;
1436   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1437   if (ts->ops->prestep) {
1438     PetscStackPush("TS PreStep function");
1439     CHKMEMQ;
1440     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1441     CHKMEMQ;
1442     PetscStackPop;
1443   }
1444   PetscFunctionReturn(0);
1445 }
1446 
1447 #undef __FUNCT__
1448 #define __FUNCT__ "TSDefaultPreStep"
1449 /*@
1450   TSDefaultPreStep - The default pre-stepping function which does nothing.
1451 
1452   Collective on TS
1453 
1454   Input Parameters:
1455 . ts  - The TS context obtained from TSCreate()
1456 
1457   Level: developer
1458 
1459 .keywords: TS, timestep
1460 @*/
1461 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts)
1462 {
1463   PetscFunctionBegin;
1464   PetscFunctionReturn(0);
1465 }
1466 
1467 #undef __FUNCT__
1468 #define __FUNCT__ "TSSetPostStep"
1469 /*@C
1470   TSSetPostStep - Sets the general-purpose function
1471   called once at the end of each time step.
1472 
1473   Logically Collective on TS
1474 
1475   Input Parameters:
1476 + ts   - The TS context obtained from TSCreate()
1477 - func - The function
1478 
1479   Calling sequence of func:
1480 . func (TS ts);
1481 
1482   Level: intermediate
1483 
1484 .keywords: TS, timestep
1485 @*/
1486 PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1487 {
1488   PetscFunctionBegin;
1489   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1490   ts->ops->poststep = func;
1491   PetscFunctionReturn(0);
1492 }
1493 
1494 #undef __FUNCT__
1495 #define __FUNCT__ "TSPostStep"
1496 /*@C
1497   TSPostStep - Runs the user-defined post-step function.
1498 
1499   Collective on TS
1500 
1501   Input Parameters:
1502 . ts   - The TS context obtained from TSCreate()
1503 
1504   Notes:
1505   TSPostStep() is typically used within time stepping implementations,
1506   so most users would not generally call this routine themselves.
1507 
1508   Level: developer
1509 
1510 .keywords: TS, timestep
1511 @*/
1512 PetscErrorCode PETSCTS_DLLEXPORT TSPostStep(TS ts)
1513 {
1514   PetscErrorCode ierr;
1515 
1516   PetscFunctionBegin;
1517   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1518   if (ts->ops->poststep) {
1519     PetscStackPush("TS PostStep function");
1520     CHKMEMQ;
1521     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1522     CHKMEMQ;
1523     PetscStackPop;
1524   }
1525   PetscFunctionReturn(0);
1526 }
1527 
1528 #undef __FUNCT__
1529 #define __FUNCT__ "TSDefaultPostStep"
1530 /*@
1531   TSDefaultPostStep - The default post-stepping function which does nothing.
1532 
1533   Collective on TS
1534 
1535   Input Parameters:
1536 . ts  - The TS context obtained from TSCreate()
1537 
1538   Level: developer
1539 
1540 .keywords: TS, timestep
1541 @*/
1542 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts)
1543 {
1544   PetscFunctionBegin;
1545   PetscFunctionReturn(0);
1546 }
1547 
1548 /* ------------ Routines to set performance monitoring options ----------- */
1549 
1550 #undef __FUNCT__
1551 #define __FUNCT__ "TSMonitorSet"
1552 /*@C
1553    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1554    timestep to display the iteration's  progress.
1555 
1556    Logically Collective on TS
1557 
1558    Input Parameters:
1559 +  ts - the TS context obtained from TSCreate()
1560 .  func - monitoring routine
1561 .  mctx - [optional] user-defined context for private data for the
1562              monitor routine (use PETSC_NULL if no context is desired)
1563 -  monitordestroy - [optional] routine that frees monitor context
1564           (may be PETSC_NULL)
1565 
1566    Calling sequence of func:
1567 $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1568 
1569 +    ts - the TS context
1570 .    steps - iteration number
1571 .    time - current time
1572 .    x - current iterate
1573 -    mctx - [optional] monitoring context
1574 
1575    Notes:
1576    This routine adds an additional monitor to the list of monitors that
1577    already has been loaded.
1578 
1579    Fortran notes: Only a single monitor function can be set for each TS object
1580 
1581    Level: intermediate
1582 
1583 .keywords: TS, timestep, set, monitor
1584 
1585 .seealso: TSMonitorDefault(), TSMonitorCancel()
1586 @*/
1587 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1588 {
1589   PetscFunctionBegin;
1590   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1591   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1592   ts->monitor[ts->numbermonitors]           = monitor;
1593   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1594   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1595   PetscFunctionReturn(0);
1596 }
1597 
1598 #undef __FUNCT__
1599 #define __FUNCT__ "TSMonitorCancel"
1600 /*@C
1601    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1602 
1603    Logically Collective on TS
1604 
1605    Input Parameters:
1606 .  ts - the TS context obtained from TSCreate()
1607 
1608    Notes:
1609    There is no way to remove a single, specific monitor.
1610 
1611    Level: intermediate
1612 
1613 .keywords: TS, timestep, set, monitor
1614 
1615 .seealso: TSMonitorDefault(), TSMonitorSet()
1616 @*/
1617 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts)
1618 {
1619   PetscErrorCode ierr;
1620   PetscInt       i;
1621 
1622   PetscFunctionBegin;
1623   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1624   for (i=0; i<ts->numbermonitors; i++) {
1625     if (ts->mdestroy[i]) {
1626       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
1627     }
1628   }
1629   ts->numbermonitors = 0;
1630   PetscFunctionReturn(0);
1631 }
1632 
1633 #undef __FUNCT__
1634 #define __FUNCT__ "TSMonitorDefault"
1635 /*@
1636    TSMonitorDefault - Sets the Default monitor
1637 
1638    Level: intermediate
1639 
1640 .keywords: TS, set, monitor
1641 
1642 .seealso: TSMonitorDefault(), TSMonitorSet()
1643 @*/
1644 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1645 {
1646   PetscErrorCode          ierr;
1647   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
1648 
1649   PetscFunctionBegin;
1650   if (!ctx) {
1651     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
1652   }
1653   ierr = PetscViewerASCIIMonitorPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1654   if (!ctx) {
1655     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
1656   }
1657   PetscFunctionReturn(0);
1658 }
1659 
1660 #undef __FUNCT__
1661 #define __FUNCT__ "TSStep"
1662 /*@
1663    TSStep - Steps the requested number of timesteps.
1664 
1665    Collective on TS
1666 
1667    Input Parameter:
1668 .  ts - the TS context obtained from TSCreate()
1669 
1670    Output Parameters:
1671 +  steps - number of iterations until termination
1672 -  ptime - time until termination
1673 
1674    Level: beginner
1675 
1676 .keywords: TS, timestep, solve
1677 
1678 .seealso: TSCreate(), TSSetUp(), TSDestroy()
1679 @*/
1680 PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1681 {
1682   PetscErrorCode ierr;
1683 
1684   PetscFunctionBegin;
1685   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1686   if (!ts->setupcalled) {
1687     ierr = TSSetUp(ts);CHKERRQ(ierr);
1688   }
1689 
1690   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1691   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1692   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1693 
1694   if (!PetscPreLoadingOn) {
1695     ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr);
1696   }
1697   PetscFunctionReturn(0);
1698 }
1699 
1700 #undef __FUNCT__
1701 #define __FUNCT__ "TSSolve"
1702 /*@
1703    TSSolve - Steps the requested number of timesteps.
1704 
1705    Collective on TS
1706 
1707    Input Parameter:
1708 +  ts - the TS context obtained from TSCreate()
1709 -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
1710 
1711    Level: beginner
1712 
1713 .keywords: TS, timestep, solve
1714 
1715 .seealso: TSCreate(), TSSetSolution(), TSStep()
1716 @*/
1717 PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x)
1718 {
1719   PetscInt       steps;
1720   PetscReal      ptime;
1721   PetscErrorCode ierr;
1722 
1723   PetscFunctionBegin;
1724   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1725   /* set solution vector if provided */
1726   if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); }
1727   /* reset time step and iteration counters */
1728   ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0;
1729   /* steps the requested number of timesteps. */
1730   ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr);
1731   PetscFunctionReturn(0);
1732 }
1733 
1734 #undef __FUNCT__
1735 #define __FUNCT__ "TSMonitor"
1736 /*
1737      Runs the user provided monitor routines, if they exists.
1738 */
1739 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1740 {
1741   PetscErrorCode ierr;
1742   PetscInt       i,n = ts->numbermonitors;
1743 
1744   PetscFunctionBegin;
1745   for (i=0; i<n; i++) {
1746     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1747   }
1748   PetscFunctionReturn(0);
1749 }
1750 
1751 /* ------------------------------------------------------------------------*/
1752 
1753 #undef __FUNCT__
1754 #define __FUNCT__ "TSMonitorLGCreate"
1755 /*@C
1756    TSMonitorLGCreate - Creates a line graph context for use with
1757    TS to monitor convergence of preconditioned residual norms.
1758 
1759    Collective on TS
1760 
1761    Input Parameters:
1762 +  host - the X display to open, or null for the local machine
1763 .  label - the title to put in the title bar
1764 .  x, y - the screen coordinates of the upper left coordinate of the window
1765 -  m, n - the screen width and height in pixels
1766 
1767    Output Parameter:
1768 .  draw - the drawing context
1769 
1770    Options Database Key:
1771 .  -ts_monitor_draw - automatically sets line graph monitor
1772 
1773    Notes:
1774    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1775 
1776    Level: intermediate
1777 
1778 .keywords: TS, monitor, line graph, residual, seealso
1779 
1780 .seealso: TSMonitorLGDestroy(), TSMonitorSet()
1781 
1782 @*/
1783 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1784 {
1785   PetscDraw      win;
1786   PetscErrorCode ierr;
1787 
1788   PetscFunctionBegin;
1789   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1790   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1791   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1792   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1793 
1794   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1795   PetscFunctionReturn(0);
1796 }
1797 
1798 #undef __FUNCT__
1799 #define __FUNCT__ "TSMonitorLG"
1800 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1801 {
1802   PetscDrawLG    lg = (PetscDrawLG) monctx;
1803   PetscReal      x,y = ptime;
1804   PetscErrorCode ierr;
1805 
1806   PetscFunctionBegin;
1807   if (!monctx) {
1808     MPI_Comm    comm;
1809     PetscViewer viewer;
1810 
1811     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1812     viewer = PETSC_VIEWER_DRAW_(comm);
1813     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
1814   }
1815 
1816   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1817   x = (PetscReal)n;
1818   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1819   if (n < 20 || (n % 5)) {
1820     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1821   }
1822   PetscFunctionReturn(0);
1823 }
1824 
1825 #undef __FUNCT__
1826 #define __FUNCT__ "TSMonitorLGDestroy"
1827 /*@C
1828    TSMonitorLGDestroy - Destroys a line graph context that was created
1829    with TSMonitorLGCreate().
1830 
1831    Collective on PetscDrawLG
1832 
1833    Input Parameter:
1834 .  draw - the drawing context
1835 
1836    Level: intermediate
1837 
1838 .keywords: TS, monitor, line graph, destroy
1839 
1840 .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1841 @*/
1842 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg)
1843 {
1844   PetscDraw      draw;
1845   PetscErrorCode ierr;
1846 
1847   PetscFunctionBegin;
1848   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1849   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1850   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1851   PetscFunctionReturn(0);
1852 }
1853 
1854 #undef __FUNCT__
1855 #define __FUNCT__ "TSGetTime"
1856 /*@
1857    TSGetTime - Gets the current time.
1858 
1859    Not Collective
1860 
1861    Input Parameter:
1862 .  ts - the TS context obtained from TSCreate()
1863 
1864    Output Parameter:
1865 .  t  - the current time
1866 
1867    Level: beginner
1868 
1869 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1870 
1871 .keywords: TS, get, time
1872 @*/
1873 PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t)
1874 {
1875   PetscFunctionBegin;
1876   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1877   PetscValidDoublePointer(t,2);
1878   *t = ts->ptime;
1879   PetscFunctionReturn(0);
1880 }
1881 
1882 #undef __FUNCT__
1883 #define __FUNCT__ "TSSetTime"
1884 /*@
1885    TSSetTime - Allows one to reset the time.
1886 
1887    Logically Collective on TS
1888 
1889    Input Parameters:
1890 +  ts - the TS context obtained from TSCreate()
1891 -  time - the time
1892 
1893    Level: intermediate
1894 
1895 .seealso: TSGetTime(), TSSetDuration()
1896 
1897 .keywords: TS, set, time
1898 @*/
1899 PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t)
1900 {
1901   PetscFunctionBegin;
1902   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1903   PetscValidLogicalCollectiveReal(ts,t,2);
1904   ts->ptime = t;
1905   PetscFunctionReturn(0);
1906 }
1907 
1908 #undef __FUNCT__
1909 #define __FUNCT__ "TSSetOptionsPrefix"
1910 /*@C
1911    TSSetOptionsPrefix - Sets the prefix used for searching for all
1912    TS options in the database.
1913 
1914    Logically Collective on TS
1915 
1916    Input Parameter:
1917 +  ts     - The TS context
1918 -  prefix - The prefix to prepend to all option names
1919 
1920    Notes:
1921    A hyphen (-) must NOT be given at the beginning of the prefix name.
1922    The first character of all runtime options is AUTOMATICALLY the
1923    hyphen.
1924 
1925    Level: advanced
1926 
1927 .keywords: TS, set, options, prefix, database
1928 
1929 .seealso: TSSetFromOptions()
1930 
1931 @*/
1932 PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[])
1933 {
1934   PetscErrorCode ierr;
1935 
1936   PetscFunctionBegin;
1937   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1938   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1939   switch(ts->problem_type) {
1940     case TS_NONLINEAR:
1941       if (ts->snes) {
1942         ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1943       }
1944       break;
1945     case TS_LINEAR:
1946       if (ts->ksp) {
1947         ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1948       }
1949       break;
1950   }
1951   PetscFunctionReturn(0);
1952 }
1953 
1954 
1955 #undef __FUNCT__
1956 #define __FUNCT__ "TSAppendOptionsPrefix"
1957 /*@C
1958    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1959    TS options in the database.
1960 
1961    Logically Collective on TS
1962 
1963    Input Parameter:
1964 +  ts     - The TS context
1965 -  prefix - The prefix to prepend to all option names
1966 
1967    Notes:
1968    A hyphen (-) must NOT be given at the beginning of the prefix name.
1969    The first character of all runtime options is AUTOMATICALLY the
1970    hyphen.
1971 
1972    Level: advanced
1973 
1974 .keywords: TS, append, options, prefix, database
1975 
1976 .seealso: TSGetOptionsPrefix()
1977 
1978 @*/
1979 PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[])
1980 {
1981   PetscErrorCode ierr;
1982 
1983   PetscFunctionBegin;
1984   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1985   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1986   switch(ts->problem_type) {
1987     case TS_NONLINEAR:
1988       if (ts->snes) {
1989         ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1990       }
1991       break;
1992     case TS_LINEAR:
1993       if (ts->ksp) {
1994         ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1995       }
1996       break;
1997   }
1998   PetscFunctionReturn(0);
1999 }
2000 
2001 #undef __FUNCT__
2002 #define __FUNCT__ "TSGetOptionsPrefix"
2003 /*@C
2004    TSGetOptionsPrefix - Sets the prefix used for searching for all
2005    TS options in the database.
2006 
2007    Not Collective
2008 
2009    Input Parameter:
2010 .  ts - The TS context
2011 
2012    Output Parameter:
2013 .  prefix - A pointer to the prefix string used
2014 
2015    Notes: On the fortran side, the user should pass in a string 'prifix' of
2016    sufficient length to hold the prefix.
2017 
2018    Level: intermediate
2019 
2020 .keywords: TS, get, options, prefix, database
2021 
2022 .seealso: TSAppendOptionsPrefix()
2023 @*/
2024 PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[])
2025 {
2026   PetscErrorCode ierr;
2027 
2028   PetscFunctionBegin;
2029   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2030   PetscValidPointer(prefix,2);
2031   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2032   PetscFunctionReturn(0);
2033 }
2034 
2035 #undef __FUNCT__
2036 #define __FUNCT__ "TSGetRHSJacobian"
2037 /*@C
2038    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2039 
2040    Not Collective, but parallel objects are returned if TS is parallel
2041 
2042    Input Parameter:
2043 .  ts  - The TS context obtained from TSCreate()
2044 
2045    Output Parameters:
2046 +  J   - The Jacobian J of F, where U_t = F(U,t)
2047 .  M   - The preconditioner matrix, usually the same as J
2048 - ctx - User-defined context for Jacobian evaluation routine
2049 
2050    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2051 
2052    Level: intermediate
2053 
2054 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2055 
2056 .keywords: TS, timestep, get, matrix, Jacobian
2057 @*/
2058 PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
2059 {
2060   PetscFunctionBegin;
2061   if (J) *J = ts->Arhs;
2062   if (M) *M = ts->B;
2063   if (ctx) *ctx = ts->jacP;
2064   PetscFunctionReturn(0);
2065 }
2066 
2067 #undef __FUNCT__
2068 #define __FUNCT__ "TSGetIJacobian"
2069 /*@C
2070    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
2071 
2072    Not Collective, but parallel objects are returned if TS is parallel
2073 
2074    Input Parameter:
2075 .  ts  - The TS context obtained from TSCreate()
2076 
2077    Output Parameters:
2078 +  A   - The Jacobian of F(t,U,U_t)
2079 .  B   - The preconditioner matrix, often the same as A
2080 .  f   - The function to compute the matrices
2081 - ctx - User-defined context for Jacobian evaluation routine
2082 
2083    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2084 
2085    Level: advanced
2086 
2087 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2088 
2089 .keywords: TS, timestep, get, matrix, Jacobian
2090 @*/
2091 PetscErrorCode PETSCTS_DLLEXPORT TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2092 {
2093   PetscFunctionBegin;
2094   if (A) *A = ts->A;
2095   if (B) *B = ts->B;
2096   if (f) *f = ts->ops->ijacobian;
2097   if (ctx) *ctx = ts->jacP;
2098   PetscFunctionReturn(0);
2099 }
2100 
2101 #undef __FUNCT__
2102 #define __FUNCT__ "TSMonitorSolution"
2103 /*@C
2104    TSMonitorSolution - Monitors progress of the TS solvers by calling
2105    VecView() for the solution at each timestep
2106 
2107    Collective on TS
2108 
2109    Input Parameters:
2110 +  ts - the TS context
2111 .  step - current time-step
2112 .  ptime - current time
2113 -  dummy - either a viewer or PETSC_NULL
2114 
2115    Level: intermediate
2116 
2117 .keywords: TS,  vector, monitor, view
2118 
2119 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2120 @*/
2121 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2122 {
2123   PetscErrorCode ierr;
2124   PetscViewer    viewer = (PetscViewer) dummy;
2125 
2126   PetscFunctionBegin;
2127   if (!dummy) {
2128     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2129   }
2130   ierr = VecView(x,viewer);CHKERRQ(ierr);
2131   PetscFunctionReturn(0);
2132 }
2133 
2134 
2135 #undef __FUNCT__
2136 #define __FUNCT__ "TSSetDM"
2137 /*@
2138    TSSetDM - Sets the DM that may be used by some preconditioners
2139 
2140    Logically Collective on TS and DM
2141 
2142    Input Parameters:
2143 +  ts - the preconditioner context
2144 -  dm - the dm
2145 
2146    Level: intermediate
2147 
2148 
2149 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2150 @*/
2151 PetscErrorCode PETSCTS_DLLEXPORT TSSetDM(TS ts,DM dm)
2152 {
2153   PetscErrorCode ierr;
2154 
2155   PetscFunctionBegin;
2156   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2157   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
2158   if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);}
2159   ts->dm = dm;
2160   if (ts->snes) {
2161     ierr = SNESSetDM(ts->snes,dm);CHKERRQ(ierr);
2162   }
2163   if (ts->ksp) {
2164     ierr = KSPSetDM(ts->ksp,dm);CHKERRQ(ierr);
2165     ierr = KSPSetDMActive(ts->ksp,PETSC_FALSE);CHKERRQ(ierr);
2166   }
2167   PetscFunctionReturn(0);
2168 }
2169 
2170 #undef __FUNCT__
2171 #define __FUNCT__ "TSGetDM"
2172 /*@
2173    TSGetDM - Gets the DM that may be used by some preconditioners
2174 
2175    Not Collective
2176 
2177    Input Parameter:
2178 . ts - the preconditioner context
2179 
2180    Output Parameter:
2181 .  dm - the dm
2182 
2183    Level: intermediate
2184 
2185 
2186 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2187 @*/
2188 PetscErrorCode PETSCTS_DLLEXPORT TSGetDM(TS ts,DM *dm)
2189 {
2190   PetscFunctionBegin;
2191   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2192   *dm = ts->dm;
2193   PetscFunctionReturn(0);
2194 }
2195 
2196 #undef __FUNCT__
2197 #define __FUNCT__ "SNESTSFormFunction"
2198 /*@
2199    SNESTSFormFunction - Function to evaluate nonlinear residual
2200 
2201    Logically Collective on SNES
2202 
2203    Input Parameter:
2204 + snes - nonlinear solver
2205 . X - the current state at which to evaluate the residual
2206 - ctx - user context, must be a TS
2207 
2208    Output Parameter:
2209 . F - the nonlinear residual
2210 
2211    Notes:
2212    This function is not normally called by users and is automatically registered with the SNES used by TS.
2213    It is most frequently passed to MatFDColoringSetFunction().
2214 
2215    Level: advanced
2216 
2217 .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2218 @*/
2219 PetscErrorCode PETSCTS_DLLEXPORT SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2220 {
2221   TS ts = (TS)ctx;
2222   PetscErrorCode ierr;
2223 
2224   PetscFunctionBegin;
2225   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2226   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2227   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
2228   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
2229   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
2230   PetscFunctionReturn(0);
2231 }
2232 
2233 #undef __FUNCT__
2234 #define __FUNCT__ "SNESTSFormJacobian"
2235 /*@
2236    SNESTSFormJacobian - Function to evaluate the Jacobian
2237 
2238    Collective on SNES
2239 
2240    Input Parameter:
2241 + snes - nonlinear solver
2242 . X - the current state at which to evaluate the residual
2243 - ctx - user context, must be a TS
2244 
2245    Output Parameter:
2246 + A - the Jacobian
2247 . B - the preconditioning matrix (may be the same as A)
2248 - flag - indicates any structure change in the matrix
2249 
2250    Notes:
2251    This function is not normally called by users and is automatically registered with the SNES used by TS.
2252 
2253    Level: developer
2254 
2255 .seealso: SNESSetJacobian()
2256 @*/
2257 PetscErrorCode PETSCTS_DLLEXPORT SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2258 {
2259   TS ts = (TS)ctx;
2260   PetscErrorCode ierr;
2261 
2262   PetscFunctionBegin;
2263   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2264   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2265   PetscValidPointer(A,3);
2266   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
2267   PetscValidPointer(B,4);
2268   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
2269   PetscValidPointer(flag,5);
2270   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
2271   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
2272   PetscFunctionReturn(0);
2273 }
2274