xref: /petsc/src/ts/interface/ts.c (revision 97c4aaa0041b17e0788ef5380763993e79c96a23)
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   PetscBool      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   PetscBool               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   PetscBool      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   PetscBool      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 = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr);
855     if (ts->ops->view) {
856       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
857       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
858       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
859     }
860     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
861     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
862     if (ts->problem_type == TS_NONLINEAR) {
863       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
864     }
865     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
866   } else if (isstring) {
867     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
868     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
869   }
870   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
871   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
872   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
873   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
874   PetscFunctionReturn(0);
875 }
876 
877 
878 #undef __FUNCT__
879 #define __FUNCT__ "TSSetApplicationContext"
880 /*@C
881    TSSetApplicationContext - Sets an optional user-defined context for
882    the timesteppers.
883 
884    Logically Collective on TS
885 
886    Input Parameters:
887 +  ts - the TS context obtained from TSCreate()
888 -  usrP - optional user context
889 
890    Level: intermediate
891 
892 .keywords: TS, timestep, set, application, context
893 
894 .seealso: TSGetApplicationContext()
895 @*/
896 PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP)
897 {
898   PetscFunctionBegin;
899   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
900   ts->user = usrP;
901   PetscFunctionReturn(0);
902 }
903 
904 #undef __FUNCT__
905 #define __FUNCT__ "TSGetApplicationContext"
906 /*@C
907     TSGetApplicationContext - Gets the user-defined context for the
908     timestepper.
909 
910     Not Collective
911 
912     Input Parameter:
913 .   ts - the TS context obtained from TSCreate()
914 
915     Output Parameter:
916 .   usrP - user context
917 
918     Level: intermediate
919 
920 .keywords: TS, timestep, get, application, context
921 
922 .seealso: TSSetApplicationContext()
923 @*/
924 PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP)
925 {
926   PetscFunctionBegin;
927   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
928   *usrP = ts->user;
929   PetscFunctionReturn(0);
930 }
931 
932 #undef __FUNCT__
933 #define __FUNCT__ "TSGetTimeStepNumber"
934 /*@
935    TSGetTimeStepNumber - Gets the current number of timesteps.
936 
937    Not Collective
938 
939    Input Parameter:
940 .  ts - the TS context obtained from TSCreate()
941 
942    Output Parameter:
943 .  iter - number steps so far
944 
945    Level: intermediate
946 
947 .keywords: TS, timestep, get, iteration, number
948 @*/
949 PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter)
950 {
951   PetscFunctionBegin;
952   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
953   PetscValidIntPointer(iter,2);
954   *iter = ts->steps;
955   PetscFunctionReturn(0);
956 }
957 
958 #undef __FUNCT__
959 #define __FUNCT__ "TSSetInitialTimeStep"
960 /*@
961    TSSetInitialTimeStep - Sets the initial timestep to be used,
962    as well as the initial time.
963 
964    Logically Collective on TS
965 
966    Input Parameters:
967 +  ts - the TS context obtained from TSCreate()
968 .  initial_time - the initial time
969 -  time_step - the size of the timestep
970 
971    Level: intermediate
972 
973 .seealso: TSSetTimeStep(), TSGetTimeStep()
974 
975 .keywords: TS, set, initial, timestep
976 @*/
977 PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
978 {
979   PetscFunctionBegin;
980   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
981   ts->time_step         = time_step;
982   ts->initial_time_step = time_step;
983   ts->ptime             = initial_time;
984   PetscFunctionReturn(0);
985 }
986 
987 #undef __FUNCT__
988 #define __FUNCT__ "TSSetTimeStep"
989 /*@
990    TSSetTimeStep - Allows one to reset the timestep at any time,
991    useful for simple pseudo-timestepping codes.
992 
993    Logically Collective on TS
994 
995    Input Parameters:
996 +  ts - the TS context obtained from TSCreate()
997 -  time_step - the size of the timestep
998 
999    Level: intermediate
1000 
1001 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1002 
1003 .keywords: TS, set, timestep
1004 @*/
1005 PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step)
1006 {
1007   PetscFunctionBegin;
1008   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1009   PetscValidLogicalCollectiveReal(ts,time_step,2);
1010   ts->time_step = time_step;
1011   PetscFunctionReturn(0);
1012 }
1013 
1014 #undef __FUNCT__
1015 #define __FUNCT__ "TSGetTimeStep"
1016 /*@
1017    TSGetTimeStep - Gets the current timestep size.
1018 
1019    Not Collective
1020 
1021    Input Parameter:
1022 .  ts - the TS context obtained from TSCreate()
1023 
1024    Output Parameter:
1025 .  dt - the current timestep size
1026 
1027    Level: intermediate
1028 
1029 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1030 
1031 .keywords: TS, get, timestep
1032 @*/
1033 PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt)
1034 {
1035   PetscFunctionBegin;
1036   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1037   PetscValidDoublePointer(dt,2);
1038   *dt = ts->time_step;
1039   PetscFunctionReturn(0);
1040 }
1041 
1042 #undef __FUNCT__
1043 #define __FUNCT__ "TSGetSolution"
1044 /*@
1045    TSGetSolution - Returns the solution at the present timestep. It
1046    is valid to call this routine inside the function that you are evaluating
1047    in order to move to the new timestep. This vector not changed until
1048    the solution at the next timestep has been calculated.
1049 
1050    Not Collective, but Vec returned is parallel if TS is parallel
1051 
1052    Input Parameter:
1053 .  ts - the TS context obtained from TSCreate()
1054 
1055    Output Parameter:
1056 .  v - the vector containing the solution
1057 
1058    Level: intermediate
1059 
1060 .seealso: TSGetTimeStep()
1061 
1062 .keywords: TS, timestep, get, solution
1063 @*/
1064 PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v)
1065 {
1066   PetscFunctionBegin;
1067   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1068   PetscValidPointer(v,2);
1069   *v = ts->vec_sol_always;
1070   PetscFunctionReturn(0);
1071 }
1072 
1073 /* ----- Routines to initialize and destroy a timestepper ---- */
1074 #undef __FUNCT__
1075 #define __FUNCT__ "TSSetProblemType"
1076 /*@
1077   TSSetProblemType - Sets the type of problem to be solved.
1078 
1079   Not collective
1080 
1081   Input Parameters:
1082 + ts   - The TS
1083 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1084 .vb
1085          U_t = A U
1086          U_t = A(t) U
1087          U_t = F(t,U)
1088 .ve
1089 
1090    Level: beginner
1091 
1092 .keywords: TS, problem type
1093 .seealso: TSSetUp(), TSProblemType, TS
1094 @*/
1095 PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type)
1096 {
1097   PetscFunctionBegin;
1098   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1099   ts->problem_type = type;
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 #undef __FUNCT__
1104 #define __FUNCT__ "TSGetProblemType"
1105 /*@C
1106   TSGetProblemType - Gets the type of problem to be solved.
1107 
1108   Not collective
1109 
1110   Input Parameter:
1111 . ts   - The TS
1112 
1113   Output Parameter:
1114 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1115 .vb
1116          U_t = A U
1117          U_t = A(t) U
1118          U_t = F(t,U)
1119 .ve
1120 
1121    Level: beginner
1122 
1123 .keywords: TS, problem type
1124 .seealso: TSSetUp(), TSProblemType, TS
1125 @*/
1126 PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type)
1127 {
1128   PetscFunctionBegin;
1129   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1130   PetscValidIntPointer(type,2);
1131   *type = ts->problem_type;
1132   PetscFunctionReturn(0);
1133 }
1134 
1135 #undef __FUNCT__
1136 #define __FUNCT__ "TSSetUp"
1137 /*@
1138    TSSetUp - Sets up the internal data structures for the later use
1139    of a timestepper.
1140 
1141    Collective on TS
1142 
1143    Input Parameter:
1144 .  ts - the TS context obtained from TSCreate()
1145 
1146    Notes:
1147    For basic use of the TS solvers the user need not explicitly call
1148    TSSetUp(), since these actions will automatically occur during
1149    the call to TSStep().  However, if one wishes to control this
1150    phase separately, TSSetUp() should be called after TSCreate()
1151    and optional routines of the form TSSetXXX(), but before TSStep().
1152 
1153    Level: advanced
1154 
1155 .keywords: TS, timestep, setup
1156 
1157 .seealso: TSCreate(), TSStep(), TSDestroy()
1158 @*/
1159 PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts)
1160 {
1161   PetscErrorCode ierr;
1162 
1163   PetscFunctionBegin;
1164   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1165   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1166   if (!((PetscObject)ts)->type_name) {
1167     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1168   }
1169   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1170   ts->setupcalled = 1;
1171   PetscFunctionReturn(0);
1172 }
1173 
1174 #undef __FUNCT__
1175 #define __FUNCT__ "TSDestroy"
1176 /*@
1177    TSDestroy - Destroys the timestepper context that was created
1178    with TSCreate().
1179 
1180    Collective on TS
1181 
1182    Input Parameter:
1183 .  ts - the TS context obtained from TSCreate()
1184 
1185    Level: beginner
1186 
1187 .keywords: TS, timestepper, destroy
1188 
1189 .seealso: TSCreate(), TSSetUp(), TSSolve()
1190 @*/
1191 PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts)
1192 {
1193   PetscErrorCode ierr;
1194 
1195   PetscFunctionBegin;
1196   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1197   if (--((PetscObject)ts)->refct > 0) PetscFunctionReturn(0);
1198 
1199   /* if memory was published with AMS then destroy it */
1200   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
1201 
1202   if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);}
1203   if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);}
1204   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
1205   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
1206   if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);}
1207   ierr = TSMonitorCancel(ts);CHKERRQ(ierr);
1208   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 #undef __FUNCT__
1213 #define __FUNCT__ "TSGetSNES"
1214 /*@
1215    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1216    a TS (timestepper) context. Valid only for nonlinear problems.
1217 
1218    Not Collective, but SNES is parallel if TS is parallel
1219 
1220    Input Parameter:
1221 .  ts - the TS context obtained from TSCreate()
1222 
1223    Output Parameter:
1224 .  snes - the nonlinear solver context
1225 
1226    Notes:
1227    The user can then directly manipulate the SNES context to set various
1228    options, etc.  Likewise, the user can then extract and manipulate the
1229    KSP, KSP, and PC contexts as well.
1230 
1231    TSGetSNES() does not work for integrators that do not use SNES; in
1232    this case TSGetSNES() returns PETSC_NULL in snes.
1233 
1234    Level: beginner
1235 
1236 .keywords: timestep, get, SNES
1237 @*/
1238 PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes)
1239 {
1240   PetscFunctionBegin;
1241   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1242   PetscValidPointer(snes,2);
1243   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first");
1244   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
1245   *snes = ts->snes;
1246   PetscFunctionReturn(0);
1247 }
1248 
1249 #undef __FUNCT__
1250 #define __FUNCT__ "TSGetKSP"
1251 /*@
1252    TSGetKSP - Returns the KSP (linear solver) associated with
1253    a TS (timestepper) context.
1254 
1255    Not Collective, but KSP is parallel if TS is parallel
1256 
1257    Input Parameter:
1258 .  ts - the TS context obtained from TSCreate()
1259 
1260    Output Parameter:
1261 .  ksp - the nonlinear solver context
1262 
1263    Notes:
1264    The user can then directly manipulate the KSP context to set various
1265    options, etc.  Likewise, the user can then extract and manipulate the
1266    KSP and PC contexts as well.
1267 
1268    TSGetKSP() does not work for integrators that do not use KSP;
1269    in this case TSGetKSP() returns PETSC_NULL in ksp.
1270 
1271    Level: beginner
1272 
1273 .keywords: timestep, get, KSP
1274 @*/
1275 PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp)
1276 {
1277   PetscFunctionBegin;
1278   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1279   PetscValidPointer(ksp,2);
1280   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1281   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1282   *ksp = ts->ksp;
1283   PetscFunctionReturn(0);
1284 }
1285 
1286 /* ----------- Routines to set solver parameters ---------- */
1287 
1288 #undef __FUNCT__
1289 #define __FUNCT__ "TSGetDuration"
1290 /*@
1291    TSGetDuration - Gets the maximum number of timesteps to use and
1292    maximum time for iteration.
1293 
1294    Not Collective
1295 
1296    Input Parameters:
1297 +  ts       - the TS context obtained from TSCreate()
1298 .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1299 -  maxtime  - final time to iterate to, or PETSC_NULL
1300 
1301    Level: intermediate
1302 
1303 .keywords: TS, timestep, get, maximum, iterations, time
1304 @*/
1305 PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1306 {
1307   PetscFunctionBegin;
1308   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1309   if (maxsteps) {
1310     PetscValidIntPointer(maxsteps,2);
1311     *maxsteps = ts->max_steps;
1312   }
1313   if (maxtime ) {
1314     PetscValidScalarPointer(maxtime,3);
1315     *maxtime  = ts->max_time;
1316   }
1317   PetscFunctionReturn(0);
1318 }
1319 
1320 #undef __FUNCT__
1321 #define __FUNCT__ "TSSetDuration"
1322 /*@
1323    TSSetDuration - Sets the maximum number of timesteps to use and
1324    maximum time for iteration.
1325 
1326    Logically Collective on TS
1327 
1328    Input Parameters:
1329 +  ts - the TS context obtained from TSCreate()
1330 .  maxsteps - maximum number of iterations to use
1331 -  maxtime - final time to iterate to
1332 
1333    Options Database Keys:
1334 .  -ts_max_steps <maxsteps> - Sets maxsteps
1335 .  -ts_max_time <maxtime> - Sets maxtime
1336 
1337    Notes:
1338    The default maximum number of iterations is 5000. Default time is 5.0
1339 
1340    Level: intermediate
1341 
1342 .keywords: TS, timestep, set, maximum, iterations
1343 @*/
1344 PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1345 {
1346   PetscFunctionBegin;
1347   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1348   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1349   PetscValidLogicalCollectiveReal(ts,maxtime,2);
1350   ts->max_steps = maxsteps;
1351   ts->max_time  = maxtime;
1352   PetscFunctionReturn(0);
1353 }
1354 
1355 #undef __FUNCT__
1356 #define __FUNCT__ "TSSetSolution"
1357 /*@
1358    TSSetSolution - Sets the initial solution vector
1359    for use by the TS routines.
1360 
1361    Logically Collective on TS and Vec
1362 
1363    Input Parameters:
1364 +  ts - the TS context obtained from TSCreate()
1365 -  x - the solution vector
1366 
1367    Level: beginner
1368 
1369 .keywords: TS, timestep, set, solution, initial conditions
1370 @*/
1371 PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x)
1372 {
1373   PetscFunctionBegin;
1374   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1375   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1376   ts->vec_sol        = ts->vec_sol_always = x;
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 #undef __FUNCT__
1381 #define __FUNCT__ "TSSetPreStep"
1382 /*@C
1383   TSSetPreStep - Sets the general-purpose function
1384   called once at the beginning of each time step.
1385 
1386   Logically Collective on TS
1387 
1388   Input Parameters:
1389 + ts   - The TS context obtained from TSCreate()
1390 - func - The function
1391 
1392   Calling sequence of func:
1393 . func (TS ts);
1394 
1395   Level: intermediate
1396 
1397 .keywords: TS, timestep
1398 @*/
1399 PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1400 {
1401   PetscFunctionBegin;
1402   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1403   ts->ops->prestep = func;
1404   PetscFunctionReturn(0);
1405 }
1406 
1407 #undef __FUNCT__
1408 #define __FUNCT__ "TSPreStep"
1409 /*@C
1410   TSPreStep - Runs the user-defined pre-step function.
1411 
1412   Collective on TS
1413 
1414   Input Parameters:
1415 . ts   - The TS context obtained from TSCreate()
1416 
1417   Notes:
1418   TSPreStep() is typically used within time stepping implementations,
1419   so most users would not generally call this routine themselves.
1420 
1421   Level: developer
1422 
1423 .keywords: TS, timestep
1424 @*/
1425 PetscErrorCode PETSCTS_DLLEXPORT TSPreStep(TS ts)
1426 {
1427   PetscErrorCode ierr;
1428 
1429   PetscFunctionBegin;
1430   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1431   if (ts->ops->prestep) {
1432     PetscStackPush("TS PreStep function");
1433     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1434     PetscStackPop;
1435   }
1436   PetscFunctionReturn(0);
1437 }
1438 
1439 #undef __FUNCT__
1440 #define __FUNCT__ "TSDefaultPreStep"
1441 /*@
1442   TSDefaultPreStep - The default pre-stepping function which does nothing.
1443 
1444   Collective on TS
1445 
1446   Input Parameters:
1447 . ts  - The TS context obtained from TSCreate()
1448 
1449   Level: developer
1450 
1451 .keywords: TS, timestep
1452 @*/
1453 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts)
1454 {
1455   PetscFunctionBegin;
1456   PetscFunctionReturn(0);
1457 }
1458 
1459 #undef __FUNCT__
1460 #define __FUNCT__ "TSSetPostStep"
1461 /*@C
1462   TSSetPostStep - Sets the general-purpose function
1463   called once at the end of each time step.
1464 
1465   Logically Collective on TS
1466 
1467   Input Parameters:
1468 + ts   - The TS context obtained from TSCreate()
1469 - func - The function
1470 
1471   Calling sequence of func:
1472 . func (TS ts);
1473 
1474   Level: intermediate
1475 
1476 .keywords: TS, timestep
1477 @*/
1478 PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1479 {
1480   PetscFunctionBegin;
1481   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1482   ts->ops->poststep = func;
1483   PetscFunctionReturn(0);
1484 }
1485 
1486 #undef __FUNCT__
1487 #define __FUNCT__ "TSPostStep"
1488 /*@C
1489   TSPostStep - Runs the user-defined post-step function.
1490 
1491   Collective on TS
1492 
1493   Input Parameters:
1494 . ts   - The TS context obtained from TSCreate()
1495 
1496   Notes:
1497   TSPostStep() is typically used within time stepping implementations,
1498   so most users would not generally call this routine themselves.
1499 
1500   Level: developer
1501 
1502 .keywords: TS, timestep
1503 @*/
1504 PetscErrorCode PETSCTS_DLLEXPORT TSPostStep(TS ts)
1505 {
1506   PetscErrorCode ierr;
1507 
1508   PetscFunctionBegin;
1509   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1510   if (ts->ops->poststep) {
1511     PetscStackPush("TS PostStep function");
1512     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1513     PetscStackPop;
1514   }
1515   PetscFunctionReturn(0);
1516 }
1517 
1518 #undef __FUNCT__
1519 #define __FUNCT__ "TSDefaultPostStep"
1520 /*@
1521   TSDefaultPostStep - The default post-stepping function which does nothing.
1522 
1523   Collective on TS
1524 
1525   Input Parameters:
1526 . ts  - The TS context obtained from TSCreate()
1527 
1528   Level: developer
1529 
1530 .keywords: TS, timestep
1531 @*/
1532 PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts)
1533 {
1534   PetscFunctionBegin;
1535   PetscFunctionReturn(0);
1536 }
1537 
1538 /* ------------ Routines to set performance monitoring options ----------- */
1539 
1540 #undef __FUNCT__
1541 #define __FUNCT__ "TSMonitorSet"
1542 /*@C
1543    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1544    timestep to display the iteration's  progress.
1545 
1546    Logically Collective on TS
1547 
1548    Input Parameters:
1549 +  ts - the TS context obtained from TSCreate()
1550 .  func - monitoring routine
1551 .  mctx - [optional] user-defined context for private data for the
1552              monitor routine (use PETSC_NULL if no context is desired)
1553 -  monitordestroy - [optional] routine that frees monitor context
1554           (may be PETSC_NULL)
1555 
1556    Calling sequence of func:
1557 $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1558 
1559 +    ts - the TS context
1560 .    steps - iteration number
1561 .    time - current time
1562 .    x - current iterate
1563 -    mctx - [optional] monitoring context
1564 
1565    Notes:
1566    This routine adds an additional monitor to the list of monitors that
1567    already has been loaded.
1568 
1569    Fortran notes: Only a single monitor function can be set for each TS object
1570 
1571    Level: intermediate
1572 
1573 .keywords: TS, timestep, set, monitor
1574 
1575 .seealso: TSMonitorDefault(), TSMonitorCancel()
1576 @*/
1577 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1578 {
1579   PetscFunctionBegin;
1580   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1581   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1582   ts->monitor[ts->numbermonitors]           = monitor;
1583   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1584   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1585   PetscFunctionReturn(0);
1586 }
1587 
1588 #undef __FUNCT__
1589 #define __FUNCT__ "TSMonitorCancel"
1590 /*@C
1591    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1592 
1593    Logically Collective on TS
1594 
1595    Input Parameters:
1596 .  ts - the TS context obtained from TSCreate()
1597 
1598    Notes:
1599    There is no way to remove a single, specific monitor.
1600 
1601    Level: intermediate
1602 
1603 .keywords: TS, timestep, set, monitor
1604 
1605 .seealso: TSMonitorDefault(), TSMonitorSet()
1606 @*/
1607 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts)
1608 {
1609   PetscErrorCode ierr;
1610   PetscInt       i;
1611 
1612   PetscFunctionBegin;
1613   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1614   for (i=0; i<ts->numbermonitors; i++) {
1615     if (ts->mdestroy[i]) {
1616       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
1617     }
1618   }
1619   ts->numbermonitors = 0;
1620   PetscFunctionReturn(0);
1621 }
1622 
1623 #undef __FUNCT__
1624 #define __FUNCT__ "TSMonitorDefault"
1625 /*@
1626    TSMonitorDefault - Sets the Default monitor
1627 
1628    Level: intermediate
1629 
1630 .keywords: TS, set, monitor
1631 
1632 .seealso: TSMonitorDefault(), TSMonitorSet()
1633 @*/
1634 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1635 {
1636   PetscErrorCode          ierr;
1637   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
1638 
1639   PetscFunctionBegin;
1640   if (!ctx) {
1641     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
1642   }
1643   ierr = PetscViewerASCIIMonitorPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1644   if (!ctx) {
1645     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
1646   }
1647   PetscFunctionReturn(0);
1648 }
1649 
1650 #undef __FUNCT__
1651 #define __FUNCT__ "TSStep"
1652 /*@
1653    TSStep - Steps the requested number of timesteps.
1654 
1655    Collective on TS
1656 
1657    Input Parameter:
1658 .  ts - the TS context obtained from TSCreate()
1659 
1660    Output Parameters:
1661 +  steps - number of iterations until termination
1662 -  ptime - time until termination
1663 
1664    Level: beginner
1665 
1666 .keywords: TS, timestep, solve
1667 
1668 .seealso: TSCreate(), TSSetUp(), TSDestroy()
1669 @*/
1670 PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1671 {
1672   PetscErrorCode ierr;
1673 
1674   PetscFunctionBegin;
1675   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1676   if (!ts->setupcalled) {
1677     ierr = TSSetUp(ts);CHKERRQ(ierr);
1678   }
1679 
1680   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1681   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1682   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1683 
1684   if (!PetscPreLoadingOn) {
1685     ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr);
1686   }
1687   PetscFunctionReturn(0);
1688 }
1689 
1690 #undef __FUNCT__
1691 #define __FUNCT__ "TSSolve"
1692 /*@
1693    TSSolve - Steps the requested number of timesteps.
1694 
1695    Collective on TS
1696 
1697    Input Parameter:
1698 +  ts - the TS context obtained from TSCreate()
1699 -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
1700 
1701    Level: beginner
1702 
1703 .keywords: TS, timestep, solve
1704 
1705 .seealso: TSCreate(), TSSetSolution(), TSStep()
1706 @*/
1707 PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x)
1708 {
1709   PetscInt       steps;
1710   PetscReal      ptime;
1711   PetscErrorCode ierr;
1712 
1713   PetscFunctionBegin;
1714   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1715   /* set solution vector if provided */
1716   if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); }
1717   /* reset time step and iteration counters */
1718   ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0;
1719   /* steps the requested number of timesteps. */
1720   ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr);
1721   PetscFunctionReturn(0);
1722 }
1723 
1724 #undef __FUNCT__
1725 #define __FUNCT__ "TSMonitor"
1726 /*
1727      Runs the user provided monitor routines, if they exists.
1728 */
1729 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1730 {
1731   PetscErrorCode ierr;
1732   PetscInt       i,n = ts->numbermonitors;
1733 
1734   PetscFunctionBegin;
1735   for (i=0; i<n; i++) {
1736     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1737   }
1738   PetscFunctionReturn(0);
1739 }
1740 
1741 /* ------------------------------------------------------------------------*/
1742 
1743 #undef __FUNCT__
1744 #define __FUNCT__ "TSMonitorLGCreate"
1745 /*@C
1746    TSMonitorLGCreate - Creates a line graph context for use with
1747    TS to monitor convergence of preconditioned residual norms.
1748 
1749    Collective on TS
1750 
1751    Input Parameters:
1752 +  host - the X display to open, or null for the local machine
1753 .  label - the title to put in the title bar
1754 .  x, y - the screen coordinates of the upper left coordinate of the window
1755 -  m, n - the screen width and height in pixels
1756 
1757    Output Parameter:
1758 .  draw - the drawing context
1759 
1760    Options Database Key:
1761 .  -ts_monitor_draw - automatically sets line graph monitor
1762 
1763    Notes:
1764    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1765 
1766    Level: intermediate
1767 
1768 .keywords: TS, monitor, line graph, residual, seealso
1769 
1770 .seealso: TSMonitorLGDestroy(), TSMonitorSet()
1771 
1772 @*/
1773 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1774 {
1775   PetscDraw      win;
1776   PetscErrorCode ierr;
1777 
1778   PetscFunctionBegin;
1779   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1780   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1781   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1782   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1783 
1784   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1785   PetscFunctionReturn(0);
1786 }
1787 
1788 #undef __FUNCT__
1789 #define __FUNCT__ "TSMonitorLG"
1790 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1791 {
1792   PetscDrawLG    lg = (PetscDrawLG) monctx;
1793   PetscReal      x,y = ptime;
1794   PetscErrorCode ierr;
1795 
1796   PetscFunctionBegin;
1797   if (!monctx) {
1798     MPI_Comm    comm;
1799     PetscViewer viewer;
1800 
1801     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1802     viewer = PETSC_VIEWER_DRAW_(comm);
1803     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
1804   }
1805 
1806   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1807   x = (PetscReal)n;
1808   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1809   if (n < 20 || (n % 5)) {
1810     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1811   }
1812   PetscFunctionReturn(0);
1813 }
1814 
1815 #undef __FUNCT__
1816 #define __FUNCT__ "TSMonitorLGDestroy"
1817 /*@C
1818    TSMonitorLGDestroy - Destroys a line graph context that was created
1819    with TSMonitorLGCreate().
1820 
1821    Collective on PetscDrawLG
1822 
1823    Input Parameter:
1824 .  draw - the drawing context
1825 
1826    Level: intermediate
1827 
1828 .keywords: TS, monitor, line graph, destroy
1829 
1830 .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1831 @*/
1832 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg)
1833 {
1834   PetscDraw      draw;
1835   PetscErrorCode ierr;
1836 
1837   PetscFunctionBegin;
1838   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1839   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1840   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1841   PetscFunctionReturn(0);
1842 }
1843 
1844 #undef __FUNCT__
1845 #define __FUNCT__ "TSGetTime"
1846 /*@
1847    TSGetTime - Gets the current time.
1848 
1849    Not Collective
1850 
1851    Input Parameter:
1852 .  ts - the TS context obtained from TSCreate()
1853 
1854    Output Parameter:
1855 .  t  - the current time
1856 
1857    Level: beginner
1858 
1859 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1860 
1861 .keywords: TS, get, time
1862 @*/
1863 PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t)
1864 {
1865   PetscFunctionBegin;
1866   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1867   PetscValidDoublePointer(t,2);
1868   *t = ts->ptime;
1869   PetscFunctionReturn(0);
1870 }
1871 
1872 #undef __FUNCT__
1873 #define __FUNCT__ "TSSetTime"
1874 /*@
1875    TSSetTime - Allows one to reset the time.
1876 
1877    Logically Collective on TS
1878 
1879    Input Parameters:
1880 +  ts - the TS context obtained from TSCreate()
1881 -  time - the time
1882 
1883    Level: intermediate
1884 
1885 .seealso: TSGetTime(), TSSetDuration()
1886 
1887 .keywords: TS, set, time
1888 @*/
1889 PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t)
1890 {
1891   PetscFunctionBegin;
1892   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1893   PetscValidLogicalCollectiveReal(ts,t,2);
1894   ts->ptime = t;
1895   PetscFunctionReturn(0);
1896 }
1897 
1898 #undef __FUNCT__
1899 #define __FUNCT__ "TSSetOptionsPrefix"
1900 /*@C
1901    TSSetOptionsPrefix - Sets the prefix used for searching for all
1902    TS options in the database.
1903 
1904    Logically Collective on TS
1905 
1906    Input Parameter:
1907 +  ts     - The TS context
1908 -  prefix - The prefix to prepend to all option names
1909 
1910    Notes:
1911    A hyphen (-) must NOT be given at the beginning of the prefix name.
1912    The first character of all runtime options is AUTOMATICALLY the
1913    hyphen.
1914 
1915    Level: advanced
1916 
1917 .keywords: TS, set, options, prefix, database
1918 
1919 .seealso: TSSetFromOptions()
1920 
1921 @*/
1922 PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[])
1923 {
1924   PetscErrorCode ierr;
1925 
1926   PetscFunctionBegin;
1927   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1928   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1929   switch(ts->problem_type) {
1930     case TS_NONLINEAR:
1931       if (ts->snes) {
1932         ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1933       }
1934       break;
1935     case TS_LINEAR:
1936       if (ts->ksp) {
1937         ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1938       }
1939       break;
1940   }
1941   PetscFunctionReturn(0);
1942 }
1943 
1944 
1945 #undef __FUNCT__
1946 #define __FUNCT__ "TSAppendOptionsPrefix"
1947 /*@C
1948    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1949    TS options in the database.
1950 
1951    Logically Collective on TS
1952 
1953    Input Parameter:
1954 +  ts     - The TS context
1955 -  prefix - The prefix to prepend to all option names
1956 
1957    Notes:
1958    A hyphen (-) must NOT be given at the beginning of the prefix name.
1959    The first character of all runtime options is AUTOMATICALLY the
1960    hyphen.
1961 
1962    Level: advanced
1963 
1964 .keywords: TS, append, options, prefix, database
1965 
1966 .seealso: TSGetOptionsPrefix()
1967 
1968 @*/
1969 PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[])
1970 {
1971   PetscErrorCode ierr;
1972 
1973   PetscFunctionBegin;
1974   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1975   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1976   switch(ts->problem_type) {
1977     case TS_NONLINEAR:
1978       if (ts->snes) {
1979         ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1980       }
1981       break;
1982     case TS_LINEAR:
1983       if (ts->ksp) {
1984         ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1985       }
1986       break;
1987   }
1988   PetscFunctionReturn(0);
1989 }
1990 
1991 #undef __FUNCT__
1992 #define __FUNCT__ "TSGetOptionsPrefix"
1993 /*@C
1994    TSGetOptionsPrefix - Sets the prefix used for searching for all
1995    TS options in the database.
1996 
1997    Not Collective
1998 
1999    Input Parameter:
2000 .  ts - The TS context
2001 
2002    Output Parameter:
2003 .  prefix - A pointer to the prefix string used
2004 
2005    Notes: On the fortran side, the user should pass in a string 'prifix' of
2006    sufficient length to hold the prefix.
2007 
2008    Level: intermediate
2009 
2010 .keywords: TS, get, options, prefix, database
2011 
2012 .seealso: TSAppendOptionsPrefix()
2013 @*/
2014 PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[])
2015 {
2016   PetscErrorCode ierr;
2017 
2018   PetscFunctionBegin;
2019   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2020   PetscValidPointer(prefix,2);
2021   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2022   PetscFunctionReturn(0);
2023 }
2024 
2025 #undef __FUNCT__
2026 #define __FUNCT__ "TSGetRHSJacobian"
2027 /*@C
2028    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2029 
2030    Not Collective, but parallel objects are returned if TS is parallel
2031 
2032    Input Parameter:
2033 .  ts  - The TS context obtained from TSCreate()
2034 
2035    Output Parameters:
2036 +  J   - The Jacobian J of F, where U_t = F(U,t)
2037 .  M   - The preconditioner matrix, usually the same as J
2038 - ctx - User-defined context for Jacobian evaluation routine
2039 
2040    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2041 
2042    Level: intermediate
2043 
2044 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2045 
2046 .keywords: TS, timestep, get, matrix, Jacobian
2047 @*/
2048 PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
2049 {
2050   PetscFunctionBegin;
2051   if (J) *J = ts->Arhs;
2052   if (M) *M = ts->B;
2053   if (ctx) *ctx = ts->jacP;
2054   PetscFunctionReturn(0);
2055 }
2056 
2057 #undef __FUNCT__
2058 #define __FUNCT__ "TSGetIJacobian"
2059 /*@C
2060    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
2061 
2062    Not Collective, but parallel objects are returned if TS is parallel
2063 
2064    Input Parameter:
2065 .  ts  - The TS context obtained from TSCreate()
2066 
2067    Output Parameters:
2068 +  A   - The Jacobian of F(t,U,U_t)
2069 .  B   - The preconditioner matrix, often the same as A
2070 .  f   - The function to compute the matrices
2071 - ctx - User-defined context for Jacobian evaluation routine
2072 
2073    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2074 
2075    Level: advanced
2076 
2077 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2078 
2079 .keywords: TS, timestep, get, matrix, Jacobian
2080 @*/
2081 PetscErrorCode PETSCTS_DLLEXPORT TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2082 {
2083   PetscFunctionBegin;
2084   if (A) *A = ts->A;
2085   if (B) *B = ts->B;
2086   if (f) *f = ts->ops->ijacobian;
2087   if (ctx) *ctx = ts->jacP;
2088   PetscFunctionReturn(0);
2089 }
2090 
2091 #undef __FUNCT__
2092 #define __FUNCT__ "TSMonitorSolution"
2093 /*@C
2094    TSMonitorSolution - Monitors progress of the TS solvers by calling
2095    VecView() for the solution at each timestep
2096 
2097    Collective on TS
2098 
2099    Input Parameters:
2100 +  ts - the TS context
2101 .  step - current time-step
2102 .  ptime - current time
2103 -  dummy - either a viewer or PETSC_NULL
2104 
2105    Level: intermediate
2106 
2107 .keywords: TS,  vector, monitor, view
2108 
2109 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2110 @*/
2111 PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2112 {
2113   PetscErrorCode ierr;
2114   PetscViewer    viewer = (PetscViewer) dummy;
2115 
2116   PetscFunctionBegin;
2117   if (!dummy) {
2118     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2119   }
2120   ierr = VecView(x,viewer);CHKERRQ(ierr);
2121   PetscFunctionReturn(0);
2122 }
2123 
2124 
2125 #undef __FUNCT__
2126 #define __FUNCT__ "TSSetDM"
2127 /*@
2128    TSSetDM - Sets the DM that may be used by some preconditioners
2129 
2130    Logically Collective on TS and DM
2131 
2132    Input Parameters:
2133 +  ts - the preconditioner context
2134 -  dm - the dm
2135 
2136    Level: intermediate
2137 
2138 
2139 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2140 @*/
2141 PetscErrorCode PETSCTS_DLLEXPORT TSSetDM(TS ts,DM dm)
2142 {
2143   PetscErrorCode ierr;
2144 
2145   PetscFunctionBegin;
2146   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2147   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
2148   if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);}
2149   ts->dm = dm;
2150   if (ts->snes) {
2151     ierr = SNESSetDM(ts->snes,dm);CHKERRQ(ierr);
2152   }
2153   if (ts->ksp) {
2154     ierr = KSPSetDM(ts->ksp,dm);CHKERRQ(ierr);
2155     ierr = KSPSetDMActive(ts->ksp,PETSC_FALSE);CHKERRQ(ierr);
2156   }
2157   PetscFunctionReturn(0);
2158 }
2159 
2160 #undef __FUNCT__
2161 #define __FUNCT__ "TSGetDM"
2162 /*@
2163    TSGetDM - Gets the DM that may be used by some preconditioners
2164 
2165    Not Collective
2166 
2167    Input Parameter:
2168 . ts - the preconditioner context
2169 
2170    Output Parameter:
2171 .  dm - the dm
2172 
2173    Level: intermediate
2174 
2175 
2176 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2177 @*/
2178 PetscErrorCode PETSCTS_DLLEXPORT TSGetDM(TS ts,DM *dm)
2179 {
2180   PetscFunctionBegin;
2181   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2182   *dm = ts->dm;
2183   PetscFunctionReturn(0);
2184 }
2185 
2186 #undef __FUNCT__
2187 #define __FUNCT__ "SNESTSFormFunction"
2188 /*@
2189    SNESTSFormFunction - Function to evaluate nonlinear residual
2190 
2191    Logically Collective on SNES
2192 
2193    Input Parameter:
2194 + snes - nonlinear solver
2195 . X - the current state at which to evaluate the residual
2196 - ctx - user context, must be a TS
2197 
2198    Output Parameter:
2199 . F - the nonlinear residual
2200 
2201    Notes:
2202    This function is not normally called by users and is automatically registered with the SNES used by TS.
2203    It is most frequently passed to MatFDColoringSetFunction().
2204 
2205    Level: advanced
2206 
2207 .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2208 @*/
2209 PetscErrorCode PETSCTS_DLLEXPORT SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2210 {
2211   TS ts = (TS)ctx;
2212   PetscErrorCode ierr;
2213 
2214   PetscFunctionBegin;
2215   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2216   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2217   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
2218   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
2219   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
2220   PetscFunctionReturn(0);
2221 }
2222 
2223 #undef __FUNCT__
2224 #define __FUNCT__ "SNESTSFormJacobian"
2225 /*@
2226    SNESTSFormJacobian - Function to evaluate the Jacobian
2227 
2228    Collective on SNES
2229 
2230    Input Parameter:
2231 + snes - nonlinear solver
2232 . X - the current state at which to evaluate the residual
2233 - ctx - user context, must be a TS
2234 
2235    Output Parameter:
2236 + A - the Jacobian
2237 . B - the preconditioning matrix (may be the same as A)
2238 - flag - indicates any structure change in the matrix
2239 
2240    Notes:
2241    This function is not normally called by users and is automatically registered with the SNES used by TS.
2242 
2243    Level: developer
2244 
2245 .seealso: SNESSetJacobian()
2246 @*/
2247 PetscErrorCode PETSCTS_DLLEXPORT SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2248 {
2249   TS ts = (TS)ctx;
2250   PetscErrorCode ierr;
2251 
2252   PetscFunctionBegin;
2253   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2254   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
2255   PetscValidPointer(A,3);
2256   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
2257   PetscValidPointer(B,4);
2258   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
2259   PetscValidPointer(flag,5);
2260   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
2261   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
2262   PetscFunctionReturn(0);
2263 }
2264